CMS 3D CMS Logo

HGCalLayerClusterProducer.cc
Go to the documentation of this file.
1 #ifndef __RecoLocalCalo_HGCRecProducers_HGCalLayerClusterProducer_H__
2 #define __RecoLocalCalo_HGCRecProducers_HGCalLayerClusterProducer_H__
3 
4 // user include files
8 
15 
23 
26 
29 
32 
34 
36 public:
39  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
40 
41  void produce(edm::Event&, const edm::EventSetup&) override;
42 
43 private:
47 
49 
50  std::unique_ptr<HGCalClusteringAlgoBase> algo;
51  bool doSharing;
53 
55  double timeOffset;
56 };
57 
59 
61  : algoId(reco::CaloCluster::undefined),
62  doSharing(ps.getParameter<bool>("doSharing")),
63  detector(ps.getParameter<std::string>("detector")), // one of EE, FH, BH or "all"
64  timeClname(ps.getParameter<std::string>("timeClname")),
65  timeOffset(ps.getParameter<double>("timeOffset")) {
66  if (detector == "all") {
67  hits_ee_token = consumes<HGCRecHitCollection>(ps.getParameter<edm::InputTag>("HGCEEInput"));
68  hits_fh_token = consumes<HGCRecHitCollection>(ps.getParameter<edm::InputTag>("HGCFHInput"));
69  hits_bh_token = consumes<HGCRecHitCollection>(ps.getParameter<edm::InputTag>("HGCBHInput"));
71  } else if (detector == "EE") {
72  hits_ee_token = consumes<HGCRecHitCollection>(ps.getParameter<edm::InputTag>("HGCEEInput"));
74  } else if (detector == "FH") {
75  hits_fh_token = consumes<HGCRecHitCollection>(ps.getParameter<edm::InputTag>("HGCFHInput"));
77  } else {
78  hits_bh_token = consumes<HGCRecHitCollection>(ps.getParameter<edm::InputTag>("HGCBHInput"));
80  }
81 
82  auto pluginPSet = ps.getParameter<edm::ParameterSet>("plugin");
83  algo = HGCalLayerClusterAlgoFactory::get()->create(pluginPSet.getParameter<std::string>("type"), pluginPSet);
84  algo->setAlgoId(algoId);
85 
86  produces<std::vector<float>>("InitialLayerClustersMask");
87  produces<std::vector<reco::BasicCluster>>();
88  produces<std::vector<reco::BasicCluster>>("sharing");
89  //density
90  produces<Density>();
91  //time for layer clusters
92  produces<edm::ValueMap<float>>(timeClname);
93 }
94 
96  // hgcalLayerClusters
99  pluginDesc.addNode(edm::PluginDescription<HGCalLayerClusterAlgoFactory>("type", "CLUE", true));
100 
101  desc.add<edm::ParameterSetDescription>("plugin", pluginDesc);
102  desc.add<std::string>("detector", "all");
103  desc.add<bool>("doSharing", false);
104  desc.add<edm::InputTag>("HGCEEInput", edm::InputTag("HGCalRecHit", "HGCEERecHits"));
105  desc.add<edm::InputTag>("HGCFHInput", edm::InputTag("HGCalRecHit", "HGCHEFRecHits"));
106  desc.add<edm::InputTag>("HGCBHInput", edm::InputTag("HGCalRecHit", "HGCHEBRecHits"));
107  desc.add<std::string>("timeClname", "timeLayerCluster");
108  desc.add<double>("timeOffset", 0.0);
109  descriptions.add("hgcalLayerClusters", desc);
110 }
111 
116 
117  std::unique_ptr<std::vector<reco::BasicCluster>> clusters(new std::vector<reco::BasicCluster>),
118  clusters_sharing(new std::vector<reco::BasicCluster>);
119  auto density = std::make_unique<Density>();
120 
121  algo->reset();
122 
123  algo->getEventSetup(es);
124 
125  //make a map detid-rechit
126  // NB for the moment just host EE and FH hits
127  // timing in digi for BH not implemented for now
128  std::unordered_map<uint32_t, float> hitmap;
129 
130  switch (algoId) {
132  evt.getByToken(hits_ee_token, ee_hits);
133  algo->populate(*ee_hits);
134  for (auto const& it : *ee_hits)
135  hitmap[it.detid().rawId()] = it.time();
136  break;
138  evt.getByToken(hits_fh_token, fh_hits);
139  evt.getByToken(hits_bh_token, bh_hits);
140  if (fh_hits.isValid()) {
141  algo->populate(*fh_hits);
142  for (auto const& it : *fh_hits)
143  hitmap[it.detid().rawId()] = it.time();
144  } else if (bh_hits.isValid()) {
145  algo->populate(*bh_hits);
146  }
147  break;
149  evt.getByToken(hits_ee_token, ee_hits);
150  algo->populate(*ee_hits);
151  for (auto const& it : *ee_hits) {
152  hitmap[it.detid().rawId()] = it.time();
153  }
154  evt.getByToken(hits_fh_token, fh_hits);
155  algo->populate(*fh_hits);
156  for (auto const& it : *fh_hits) {
157  hitmap[it.detid().rawId()] = it.time();
158  }
159  evt.getByToken(hits_bh_token, bh_hits);
160  algo->populate(*bh_hits);
161  break;
162  default:
163  break;
164  }
165  algo->makeClusters();
166  *clusters = algo->getClusters(false);
167  if (doSharing)
168  *clusters_sharing = algo->getClusters(true);
169 
170  auto clusterHandle = evt.put(std::move(clusters));
171  auto clusterHandleSharing = evt.put(std::move(clusters_sharing), "sharing");
172 
173  //Keep the density
174  *density = algo->getDensity();
175  evt.put(std::move(density));
176 
177  edm::PtrVector<reco::BasicCluster> clusterPtrs, clusterPtrsSharing;
178 
179  std::vector<float> times;
180  times.reserve(clusterHandle->size());
181 
182  for (unsigned i = 0; i < clusterHandle->size(); ++i) {
183  edm::Ptr<reco::BasicCluster> ptr(clusterHandle, i);
184  clusterPtrs.push_back(ptr);
185 
186  float timeCl = -99.;
187  const reco::CaloCluster& sCl = (*clusterHandle)[i];
188  if (sCl.size() >= 3) {
189  std::vector<float> timeClhits;
190 
191  for (auto const& hit : sCl.hitsAndFractions()) {
192  auto finder = hitmap.find(hit.first);
193  if (finder == hitmap.end())
194  continue;
195 
196  //time is computed wrt 0-25ns + offset and set to -1 if no time
197  float rhTime = finder->second;
198  if (rhTime < 0.)
199  continue;
200  timeClhits.push_back(rhTime - timeOffset);
201  }
202  if (timeClhits.size() >= 3)
203  timeCl = hgcalsimclustertime::fixSizeHighestDensity(timeClhits);
204  }
205  times.push_back(timeCl);
206  }
207  std::unique_ptr<std::vector<float>> layerClustersMask(new std::vector<float>);
208  layerClustersMask->resize(clusterHandle->size(), 1.0);
209  evt.put(std::move(layerClustersMask), "InitialLayerClustersMask");
210 
211  auto timeCl = std::make_unique<edm::ValueMap<float>>();
213  filler.insert(clusterHandle, times.begin(), times.end());
214  filler.fill();
215  evt.put(std::move(timeCl), timeClname);
216 
217  if (doSharing) {
218  for (unsigned i = 0; i < clusterHandleSharing->size(); ++i) {
219  edm::Ptr<reco::BasicCluster> ptr(clusterHandleSharing, i);
220  clusterPtrsSharing.push_back(ptr);
221  }
222  }
223 }
224 
225 #endif
float fixSizeHighestDensity(std::vector< float > &t, float deltaT=0.210, float timeWidthBy=0.5)
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
void produce(edm::Event &, const edm::EventSetup &) override
reco::CaloCluster::AlgoId algoId
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
void push_back(Ptr< T > const &iPtr)
Definition: PtrVector.h:149
ParameterDescriptionNode * addNode(ParameterDescriptionNode const &node)
std::unique_ptr< HGCalClusteringAlgoBase > algo
std::map< DetId, float > Density
edm::EDGetTokenT< HGCRecHitCollection > hits_bh_token
hgcal_clustering::Density Density
Definition: HGCalCLUEAlgo.h:27
const std::vector< std::pair< DetId, float > > & hitsAndFractions() const
Definition: CaloCluster.h:209
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isValid() const
Definition: HandleBase.h:70
edm::EDGetTokenT< HGCRecHitCollection > hits_ee_token
void add(std::string const &label, ParameterSetDescription const &psetDescription)
HGCalLayerClusterProducer(const edm::ParameterSet &)
void reserve(size_type n)
Reserve space for RefVector.
Definition: PtrVectorBase.h:93
fixed size matrix
size_t size() const
size in number of hits (e.g. in crystals for ECAL)
Definition: CaloCluster.h:186
edm::EDGetTokenT< HGCRecHitCollection > hits_fh_token
def move(src, dest)
Definition: eostools.py:511