CMS 3D CMS Logo

MtdRecoClusterToSimLayerClusterAssociatorEDProducer.cc
Go to the documentation of this file.
1 // system include files
2 #include <memory>
3 #include <string>
4 
5 // user include files
15 
17 
18 //
19 // class decleration
20 //
21 
23 public:
26 
27  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
28 
29 private:
30  void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override;
31 
35 
37 };
38 
40  const edm::ParameterSet &pset) {
41  produces<reco::SimToRecoCollectionMtd>();
42  produces<reco::RecoToSimCollectionMtd>();
43 
44  btlRecoClustersToken_ = consumes<FTLClusterCollection>(pset.getParameter<edm::InputTag>("btlRecoClustersTag"));
45  etlRecoClustersToken_ = consumes<FTLClusterCollection>(pset.getParameter<edm::InputTag>("etlRecoClustersTag"));
46  simClustersToken_ = consumes<MtdSimLayerClusterCollection>(pset.getParameter<edm::InputTag>("mtdSimClustersTag"));
48  consumes<reco::MtdRecoClusterToSimLayerClusterAssociator>(pset.getParameter<edm::InputTag>("associator"));
49 }
50 
52 
53 //
54 // member functions
55 //
56 
57 // ------------ method called to produce the data ------------
60  const edm::EventSetup &iSetup) const {
61  using namespace edm;
62 
64  iEvent.getByToken(associatorToken_, theAssociator);
65 
66  edm::Handle<FTLClusterCollection> btlRecoClusters;
67  iEvent.getByToken(btlRecoClustersToken_, btlRecoClusters);
68 
69  edm::Handle<FTLClusterCollection> etlRecoClusters;
70  iEvent.getByToken(etlRecoClustersToken_, etlRecoClusters);
71 
73  iEvent.getByToken(simClustersToken_, simClusters);
74 
75  // associate reco clus to sim layer clus
76  reco::RecoToSimCollectionMtd recoToSimColl =
77  theAssociator->associateRecoToSim(btlRecoClusters, etlRecoClusters, simClusters);
78  reco::SimToRecoCollectionMtd simToRecoColl =
79  theAssociator->associateSimToReco(btlRecoClusters, etlRecoClusters, simClusters);
80 
81  auto r2s = std::make_unique<reco::RecoToSimCollectionMtd>(recoToSimColl);
82  auto s2r = std::make_unique<reco::SimToRecoCollectionMtd>(simToRecoColl);
83 
84  iEvent.put(std::move(r2s));
85  iEvent.put(std::move(s2r));
86 }
87 
90  desc.add<edm::InputTag>("associator", edm::InputTag("mtdRecoClusterToSimLayerClusterAssociatorByHits"));
91  desc.add<edm::InputTag>("mtdSimClustersTag", edm::InputTag("mix", "MergedMtdTruthLC"));
92  desc.add<edm::InputTag>("btlRecoClustersTag", edm::InputTag("mtdClusters", "FTLBarrel"));
93  desc.add<edm::InputTag>("etlRecoClustersTag", edm::InputTag("mtdClusters", "FTLEndcap"));
94 
95  cfg.add("mtdRecoClusterToSimLayerClusterAssociationDefault", desc);
96 }
97 
98 // define this as a plug-in
reco::RecoToSimCollectionMtd associateRecoToSim(const edm::Handle< FTLClusterCollection > &btlRecoClusH, const edm::Handle< FTLClusterCollection > &etlRecoClusH, const edm::Handle< MtdSimLayerClusterCollection > &simClusH) const
Associate RecoCluster to MtdSimLayerCluster.
int iEvent
Definition: GenABIO.cc:224
reco::SimToRecoCollectionMtd associateSimToReco(const edm::Handle< FTLClusterCollection > &btlRecoClusH, const edm::Handle< FTLClusterCollection > &etlRecoClusH, const edm::Handle< MtdSimLayerClusterCollection > &simClusH) const
Associate MtdSimLayerCluster to RecoCluster.
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::EDGetTokenT< reco::MtdRecoClusterToSimLayerClusterAssociator > associatorToken_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
HLT enums.
def move(src, dest)
Definition: eostools.py:511