CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
PFClusterTimeAssigner.cc
Go to the documentation of this file.
1 #ifndef __PFClusterTimeAssigner__
2 #define __PFClusterTimeAssigner__
3 
4 // user include files
7 
13 
16 
18 
20 public:
22  const edm::InputTag& clusters = conf.getParameter<edm::InputTag>("src");
23  clustersTok_ = consumes<reco::PFClusterCollection>(clusters);
24 
25  const edm::InputTag& times = conf.getParameter<edm::InputTag>("timeSrc");
26  timesTok_ = consumes<edm::ValueMap<float> >(times);
27 
28  const edm::InputTag& timeResos = conf.getParameter<edm::InputTag>("timeResoSrc");
29  timeResosTok_ = consumes<edm::ValueMap<float> >(timeResos);
30 
31  produces<reco::PFClusterCollection>();
32  }
33 
34  void produce(edm::Event& e, const edm::EventSetup& es) override;
35  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
36 
37 private:
41 };
42 
44 
46  auto clusters_out = std::make_unique<reco::PFClusterCollection>();
47 
49  e.getByToken(clustersTok_, clustersH);
50  edm::Handle<edm::ValueMap<float> > timesH, timeResosH;
51  e.getByToken(timesTok_, timesH);
52  e.getByToken(timeResosTok_, timeResosH);
53 
54  auto const& clusters = *clustersH;
55  auto const& times = *timesH;
56  auto const& timeResos = *timeResosH;
57 
58  clusters_out->reserve(clusters.size());
59  clusters_out->insert(clusters_out->end(), clusters.begin(), clusters.end());
60 
61  //build the EE->PS association
62  auto& out = *clusters_out;
63  for (unsigned i = 0; i < out.size(); ++i) {
64  edm::Ref<reco::PFClusterCollection> clusterRef(clustersH, i);
65  const float time = times[clusterRef];
66  const float timeReso = timeResos[clusterRef];
67  out[i].setTime(time, timeReso);
68  }
69 
70  e.put(std::move(clusters_out));
71 }
72 
75  desc.add<edm::InputTag>("src", edm::InputTag("particleFlowClusterECALUncorrected"));
76  desc.add<edm::InputTag>("timeSrc", edm::InputTag("ecalBarrelClusterFastTimer"));
77  desc.add<edm::InputTag>("timeResoSrc", edm::InputTag("ecalBarrelClusterFastTimer"));
78  descriptions.add("particleFlowClusterTimeAssignerDefault", desc);
79 }
80 
81 #endif
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
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_
def move
Definition: eostools.py:511
ParameterDescriptionBase * add(U const &iLabel, T const &value)
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void add(std::string const &label, ParameterSetDescription const &psetDescription)
edm::EDGetTokenT< reco::PFClusterCollection > clustersTok_