CMS 3D CMS Logo

HGCalMultiClusterProducer.cc
Go to the documentation of this file.
1 #ifndef __RecoLocalCalo_HGCRecProducers_HGCalMultiClusterProducer_H__
2 #define __RecoLocalCalo_HGCRecProducers_HGCalMultiClusterProducer_H__
3 
4 // user include files
8 
19 
23 
27 
29 
31 public:
34  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
35 
36  void produce(edm::Event&, const edm::EventSetup&) override;
37 
38 private:
44  std::unique_ptr<HGCal3DClustering> multicluster_algo;
45  bool doSharing;
47 };
48 
50 
52  : doSharing(ps.getParameter<bool>("doSharing")),
53  verbosity((HGCalClusteringAlgoBase::VerbosityLevel)ps.getUntrackedParameter<unsigned int>("verbosity", 3)) {
54  std::vector<double> multicluster_radii = ps.getParameter<std::vector<double>>("multiclusterRadii");
55  double minClusters = ps.getParameter<unsigned>("minClusters");
56  clusters_token = consumes<std::vector<reco::BasicCluster>>(ps.getParameter<edm::InputTag>("HGCLayerClusters"));
58  consumes<std::vector<reco::BasicCluster>>(ps.getParameter<edm::InputTag>("HGCLayerClustersSharing"));
59  hits_ee_token = consumes<HGCRecHitCollection>(ps.getParameter<edm::InputTag>("HGCEEInput"));
60  hits_fh_token = consumes<HGCRecHitCollection>(ps.getParameter<edm::InputTag>("HGCFHInput"));
61  hits_bh_token = consumes<HGCRecHitCollection>(ps.getParameter<edm::InputTag>("HGCBHInput"));
62  auto sumes = consumesCollector();
63 
64  multicluster_algo = std::make_unique<HGCal3DClustering>(ps, sumes, multicluster_radii, minClusters);
65 
66  produces<std::vector<reco::HGCalMultiCluster>>();
67  produces<std::vector<reco::HGCalMultiCluster>>("sharing");
68 }
69 
71  // hgcalMultiClusters
73  desc.add<edm::InputTag>("HGCLayerClusters", edm::InputTag("hgcalLayerClusters"));
74  desc.addUntracked<unsigned int>("verbosity", 3);
75  desc.add<bool>("doSharing", false);
76  desc.add<edm::InputTag>("HGCEEInput", edm::InputTag("HGCalRecHit", "HGCEERecHits"));
77  desc.add<edm::InputTag>("HGCFHInput", edm::InputTag("HGCalRecHit", "HGCHEFRecHits"));
78  desc.add<std::vector<double>>("multiclusterRadii",
79  {
80  2.0,
81  5.0,
82  5.0,
83  });
84  desc.add<edm::InputTag>("HGCBHInput", edm::InputTag("HGCalRecHit", "HGCHEBRecHits"));
85  desc.add<edm::InputTag>("HGCLayerClustersSharing", edm::InputTag("hgcalLayerClusters", "sharing"));
86  desc.add<unsigned int>("minClusters", 3);
87  descriptions.add("hgcalMultiClusters", desc);
88 }
89 
92  edm::Handle<std::vector<reco::BasicCluster>> clusterSharingHandle;
93 
94  evt.getByToken(clusters_token, clusterHandle);
95  if (doSharing)
96  evt.getByToken(clusters_sharing_token, clusterSharingHandle);
97 
98  edm::PtrVector<reco::BasicCluster> clusterPtrs, clusterPtrsSharing;
99  for (unsigned i = 0; i < clusterHandle->size(); ++i) {
100  edm::Ptr<reco::BasicCluster> ptr(clusterHandle, i);
101  clusterPtrs.push_back(ptr);
102  }
103 
104  if (doSharing) {
105  for (unsigned i = 0; i < clusterSharingHandle->size(); ++i) {
106  edm::Ptr<reco::BasicCluster> ptr(clusterSharingHandle, i);
107  clusterPtrsSharing.push_back(ptr);
108  }
109  }
110 
111  auto multiclusters = std::make_unique<std::vector<reco::HGCalMultiCluster>>();
112  auto multiclusters_sharing = std::make_unique<std::vector<reco::HGCalMultiCluster>>();
113 
114  multicluster_algo->getEvent(evt);
115  multicluster_algo->getEventSetup(es);
116 
117  *multiclusters = multicluster_algo->makeClusters(clusterPtrs);
118  if (doSharing)
119  *multiclusters_sharing = multicluster_algo->makeClusters(clusterPtrsSharing);
120  evt.put(std::move(multiclusters));
121  if (doSharing)
122  evt.put(std::move(multiclusters_sharing), "sharing");
123 }
124 
125 #endif
edm::EDGetTokenT< std::vector< reco::BasicCluster > > clusters_token
edm::EDGetTokenT< HGCRecHitCollection > hits_fh_token
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
edm::EDGetTokenT< std::vector< reco::BasicCluster > > clusters_sharing_token
edm::EDGetTokenT< HGCRecHitCollection > hits_bh_token
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void push_back(Ptr< T > const &iPtr)
Definition: PtrVector.h:149
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
HGCalClusteringAlgoBase::VerbosityLevel verbosity
edm::EDGetTokenT< HGCRecHitCollection > hits_ee_token
void produce(edm::Event &, const edm::EventSetup &) override
std::unique_ptr< HGCal3DClustering > multicluster_algo
void add(std::string const &label, ParameterSetDescription const &psetDescription)
HGCalMultiClusterProducer(const edm::ParameterSet &)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
def move(src, dest)
Definition: eostools.py:511