CMS 3D CMS Logo

TrackTimingImporter.cc
Go to the documentation of this file.
12 
13 // this doesn't actually import anything,
14 // but rather applies time stamps to tracks after they are all inserted
15 
17 public:
19  edm::ConsumesCollector& sumes) :
20  BlockElementImporterBase(conf,sumes),
21  srcTime_( sumes.consumes<edm::ValueMap<float> >(conf.getParameter<edm::InputTag>("timeValueMap")) ),
22  srcTimeError_( sumes.consumes<edm::ValueMap<float> >(conf.getParameter<edm::InputTag>("timeErrorMap")) ),
23  srcTimeGsf_( sumes.consumes<edm::ValueMap<float> >(conf.getParameter<edm::InputTag>("timeValueMapGsf")) ),
24  srcTimeErrorGsf_( sumes.consumes<edm::ValueMap<float> >(conf.getParameter<edm::InputTag>("timeErrorMapGsf")) ),
25  debug_(conf.getUntrackedParameter<bool>("debug",false)) {
26  }
27 
28  void importToBlock( const edm::Event& ,
29  ElementList& ) const override;
30 
31 private:
32 
34  const bool debug_;
35 
36 };
37 
40  "TrackTimingImporter");
41 
46 
47  edm::Handle<edm::ValueMap<float> > timeH, timeErrH, timeGsfH, timeErrGsfH;
48 
49  e.getByToken(srcTime_, timeH);
50  e.getByToken(srcTimeError_, timeErrH);
51  e.getByToken(srcTimeGsf_, timeGsfH);
52  e.getByToken(srcTimeErrorGsf_, timeErrGsfH);
53 
54  for( auto& elem : elems ) {
55  if( reco::PFBlockElement::TRACK == elem->type() ) {
56  const auto& ref = elem->trackRef();
57  if (timeH->contains(ref.id())) {
58  elem->setTime( (*timeH)[ref], (*timeErrH)[ref] );
59  }
60  if( debug_ ) {
61  edm::LogInfo("TrackTimingImporter")
62  << "Track with pT / eta " << ref->pt() << " / " << ref->eta()
63  << " has time: " << elem->time() << " +/- " << elem->timeError() << std::endl;
64  }
65  } else if ( reco::PFBlockElement::GSF == elem->type() ) {
66  const auto& ref = static_cast<const reco::PFBlockElementGsfTrack*>(elem.get())->GsftrackRef();
67  if (timeGsfH->contains(ref.id())) {
68  elem->setTime( (*timeGsfH)[ref], (*timeErrGsfH)[ref] );
69  }
70  if( debug_ ) {
71  edm::LogInfo("TrackTimingImporter")
72  << "Track with pT / eta " << ref->pt() << " / " << ref->eta()
73  << " has time: " << elem->time() << " +/- " << elem->timeError() << std::endl;
74  }
75  }
76  }
77 }
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:508
edm::EDGetTokenT< edm::ValueMap< float > > srcTimeErrorGsf_
edm::EDGetTokenT< edm::ValueMap< float > > srcTimeError_
edm::EDGetTokenT< edm::ValueMap< float > > srcTime_
def elem(elemtype, innerHTML='', html_class='', kwargs)
Definition: HTMLExport.py:18
HLT enums.
#define DEFINE_EDM_PLUGIN(factory, type, name)
std::vector< std::unique_ptr< reco::PFBlockElement > > ElementList
TrackTimingImporter(const edm::ParameterSet &conf, edm::ConsumesCollector &sumes)
void importToBlock(const edm::Event &, ElementList &) const override
edm::EDGetTokenT< edm::ValueMap< float > > srcTimeGsf_