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 
20 
24 
28 
30 
32  public:
35  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
36 
37  void produce(edm::Event&, const edm::EventSetup&) override;
38 
39  private:
45  std::unique_ptr<HGCal3DClustering> multicluster_algo;
46  bool doSharing;
48 };
49 
51 
53  doSharing(ps.getParameter<bool>("doSharing")),
54  verbosity((HGCalImagingAlgo::VerbosityLevel)ps.getUntrackedParameter<unsigned int>("verbosity",3)){
55  std::vector<double> multicluster_radii = ps.getParameter<std::vector<double> >("multiclusterRadii");
56  double minClusters = ps.getParameter<unsigned>("minClusters");
57  clusters_token = consumes<std::vector<reco::BasicCluster> >(ps.getParameter<edm::InputTag>("HGCLayerClusters"));
58  clusters_sharing_token = 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 
70 
72  // hgcalMultiClusters
74  desc.add<edm::InputTag>("HGCLayerClusters", edm::InputTag("hgcalLayerClusters"));
75  desc.addUntracked<unsigned int>("verbosity", 3);
76  desc.add<bool>("doSharing", false);
77  desc.add<edm::InputTag>("HGCEEInput", edm::InputTag("HGCalRecHit","HGCEERecHits"));
78  desc.add<edm::InputTag>("HGCFHInput", edm::InputTag("HGCalRecHit","HGCHEFRecHits"));
79  desc.add<std::vector<double>>("multiclusterRadii", {
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 
91 
93  edm::Handle<std::vector<reco::BasicCluster> > clusterSharingHandle;
94 
95 
96  evt.getByToken(clusters_token, clusterHandle);
97  if(doSharing) evt.getByToken(clusters_sharing_token, clusterSharingHandle);
98 
99 
100  edm::PtrVector<reco::BasicCluster> clusterPtrs, clusterPtrsSharing;
101  for( unsigned i = 0; i < clusterHandle->size(); ++i ) {
102  edm::Ptr<reco::BasicCluster> ptr(clusterHandle,i);
103  clusterPtrs.push_back(ptr);
104  }
105 
106  if(doSharing){
107  for( unsigned i = 0; i < clusterSharingHandle->size(); ++i ) {
108  edm::Ptr<reco::BasicCluster> ptr(clusterSharingHandle,i);
109  clusterPtrsSharing.push_back(ptr);
110  }
111  }
112 
113  auto multiclusters = std::make_unique<std::vector<reco::HGCalMultiCluster>>();
114  auto multiclusters_sharing = std::make_unique<std::vector<reco::HGCalMultiCluster>>();
115 
116 
117  multicluster_algo->getEvent(evt);
118  multicluster_algo->getEventSetup(es);
119 
120  *multiclusters = multicluster_algo->makeClusters(clusterPtrs);
121  if(doSharing)
122  *multiclusters_sharing = multicluster_algo->makeClusters(clusterPtrsSharing);
123  evt.put(std::move(multiclusters));
124  if(doSharing)
125  evt.put(std::move(multiclusters_sharing),"sharing");
126 
127 }
128 
129 #endif
T getParameter(std::string const &) const
edm::EDGetTokenT< std::vector< reco::BasicCluster > > clusters_token
edm::EDGetTokenT< HGCRecHitCollection > hits_fh_token
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:137
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
edm::EDGetTokenT< std::vector< reco::BasicCluster > > clusters_sharing_token
edm::EDGetTokenT< HGCRecHitCollection > hits_bh_token
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void push_back(Ptr< T > const &iPtr)
Definition: PtrVector.h:140
edm::EDGetTokenT< HGCRecHitCollection > hits_ee_token
void produce(edm::Event &, const edm::EventSetup &) override
HGCalImagingAlgo::VerbosityLevel verbosity
ParameterDescriptionBase * add(U const &iLabel, T const &value)
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:510