CMS 3D CMS Logo

PFClusterTimeAssigner.cc
Go to the documentation of this file.
1 #ifndef __PFClusterTimeAssigner__
2 #define __PFClusterTimeAssigner__
3 
4 // user include files
7 
13 
16 
18 
19 
21 public:
23  const edm::InputTag& clusters =
24  conf.getParameter<edm::InputTag>("src");
25  clustersTok_ = consumes<reco::PFClusterCollection>( clusters );
26 
27  const edm::InputTag& times =
28  conf.getParameter<edm::InputTag>("timeSrc");
29  timesTok_ = consumes<edm::ValueMap<float> >( times );
30 
31  const edm::InputTag& timeResos =
32  conf.getParameter<edm::InputTag>("timeResoSrc");
33  timeResosTok_ = consumes<edm::ValueMap<float> >( timeResos );
34 
35  produces<reco::PFClusterCollection>();
36  }
37 
38  void produce(edm::Event& e, const edm::EventSetup& es) override;
39  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
40 
41 private:
45 };
46 
48 
51  auto clusters_out = std::make_unique<reco::PFClusterCollection>();
52 
54  e.getByToken(clustersTok_,clustersH);
55  edm::Handle<edm::ValueMap<float> > timesH, timeResosH;
56  e.getByToken(timesTok_,timesH);
57  e.getByToken(timeResosTok_,timeResosH);
58 
59  auto const & clusters = *clustersH;
60  auto const & times = *timesH;
61  auto const & timeResos = *timeResosH;
62 
63  clusters_out->reserve(clusters.size());
64  clusters_out->insert(clusters_out->end(),
65  clusters.begin(),clusters.end());
66 
67  //build the EE->PS association
68  auto& out = *clusters_out;
69  for( unsigned i = 0; i < out.size(); ++i ) {
70 
71  edm::Ref<reco::PFClusterCollection> clusterRef(clustersH,i);
72  const float time = times[clusterRef];
73  const float timeReso = timeResos[clusterRef];
74  out[i].setTime(time, timeReso);
75 
76  }
77 
78  e.put(std::move(clusters_out));
79 }
80 
83  desc.add<edm::InputTag>("src",edm::InputTag("particleFlowClusterECALUncorrected"));
84  desc.add<edm::InputTag>("timeSrc",edm::InputTag("ecalBarrelClusterFastTimer"));
85  desc.add<edm::InputTag>("timeResoSrc",edm::InputTag("ecalBarrelClusterFastTimer"));
86  descriptions.add("particleFlowClusterTimeAssignerDefault",desc);
87 }
88 
89 #endif
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
PFClusterTimeAssigner(const edm::ParameterSet &conf)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::EDGetTokenT< edm::ValueMap< float > > timesTok_
void produce(edm::Event &e, const edm::EventSetup &es) override
edm::EDGetTokenT< edm::ValueMap< float > > timeResosTok_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
edm::EDGetTokenT< reco::PFClusterCollection > clustersTok_
def move(src, dest)
Definition: eostools.py:510