CMS 3D CMS Logo

HGCalMultiClusterProducer.cc
Go to the documentation of this file.
1 // user include files
5 
15 
19 
21 
23 public:
26  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
27 
28  void produce(edm::Event&, const edm::EventSetup&) override;
29 
30 private:
36  std::unique_ptr<HGCal3DClustering> multicluster_algo;
37  bool doSharing;
39 };
40 
42  : doSharing(ps.getParameter<bool>("doSharing")),
43  verbosity((HGCalClusteringAlgoBase::VerbosityLevel)ps.getUntrackedParameter<unsigned int>("verbosity", 3)) {
44  std::vector<double> multicluster_radii = ps.getParameter<std::vector<double>>("multiclusterRadii");
45  double minClusters = ps.getParameter<unsigned>("minClusters");
46  clusters_token = consumes<std::vector<reco::BasicCluster>>(ps.getParameter<edm::InputTag>("HGCLayerClusters"));
48  consumes<std::vector<reco::BasicCluster>>(ps.getParameter<edm::InputTag>("HGCLayerClustersSharing"));
49  hits_ee_token = consumes<HGCRecHitCollection>(ps.getParameter<edm::InputTag>("HGCEEInput"));
50  hits_fh_token = consumes<HGCRecHitCollection>(ps.getParameter<edm::InputTag>("HGCFHInput"));
51  hits_bh_token = consumes<HGCRecHitCollection>(ps.getParameter<edm::InputTag>("HGCBHInput"));
52  auto sumes = consumesCollector();
53 
54  multicluster_algo = std::make_unique<HGCal3DClustering>(ps, sumes, multicluster_radii, minClusters);
55 
56  produces<std::vector<reco::HGCalMultiCluster>>();
57  produces<std::vector<reco::HGCalMultiCluster>>("sharing");
58 }
59 
61  // hgcalMultiClusters
63  desc.add<edm::InputTag>("HGCLayerClusters", edm::InputTag("hgcalMergeLayerClusters"));
64  desc.addUntracked<unsigned int>("verbosity", 3);
65  desc.add<bool>("doSharing", false);
66  desc.add<edm::InputTag>("HGCEEInput", edm::InputTag("HGCalRecHit", "HGCEERecHits"));
67  desc.add<edm::InputTag>("HGCFHInput", edm::InputTag("HGCalRecHit", "HGCHEFRecHits"));
68  desc.add<std::vector<double>>("multiclusterRadii",
69  {
70  2.0,
71  5.0,
72  5.0,
73  });
74  desc.add<edm::InputTag>("HGCBHInput", edm::InputTag("HGCalRecHit", "HGCHEBRecHits"));
75  desc.add<edm::InputTag>("HGCLayerClustersSharing", edm::InputTag("hgcalMergeLayerClusters", "sharing"));
76  desc.add<unsigned int>("minClusters", 3);
77  descriptions.add("hgcalMultiClusters", desc);
78 }
79 
82  edm::Handle<std::vector<reco::BasicCluster>> clusterSharingHandle;
83 
84  evt.getByToken(clusters_token, clusterHandle);
85  if (doSharing)
86  evt.getByToken(clusters_sharing_token, clusterSharingHandle);
87 
88  edm::PtrVector<reco::BasicCluster> clusterPtrs, clusterPtrsSharing;
89  for (unsigned i = 0; i < clusterHandle->size(); ++i) {
90  edm::Ptr<reco::BasicCluster> ptr(clusterHandle, i);
91  clusterPtrs.push_back(ptr);
92  }
93 
94  if (doSharing) {
95  for (unsigned i = 0; i < clusterSharingHandle->size(); ++i) {
96  edm::Ptr<reco::BasicCluster> ptr(clusterSharingHandle, i);
97  clusterPtrsSharing.push_back(ptr);
98  }
99  }
100 
101  auto multiclusters = std::make_unique<std::vector<reco::HGCalMultiCluster>>();
102  auto multiclusters_sharing = std::make_unique<std::vector<reco::HGCalMultiCluster>>();
103 
104  multicluster_algo->getEvent(evt);
105  multicluster_algo->getEventSetup(es);
106 
107  *multiclusters = multicluster_algo->makeClusters(clusterPtrs);
108  if (doSharing)
109  *multiclusters_sharing = multicluster_algo->makeClusters(clusterPtrsSharing);
110  evt.put(std::move(multiclusters));
111  if (doSharing)
112  evt.put(std::move(multiclusters_sharing), "sharing");
113 }
114 
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
void push_back(Ptr< T > const &iPtr)
Definition: PtrVector.h:152
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:528
HGCalClusteringAlgoBase::VerbosityLevel verbosity
edm::EDGetTokenT< HGCRecHitCollection > hits_ee_token
void produce(edm::Event &, const edm::EventSetup &) override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const int verbosity
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)
VerbosityLevel
Definition: Common.h:71
def move(src, dest)
Definition: eostools.py:511