CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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:
48 
50 
51  std::unique_ptr<HGCalClusteringAlgoBase> algo;
52  bool doSharing;
54 
56  double timeOffset;
57  unsigned int nHitsTime;
58 };
59 
61 
63  : algoId(reco::CaloCluster::undefined),
64  doSharing(ps.getParameter<bool>("doSharing")),
65  detector(ps.getParameter<std::string>("detector")), // one of EE, FH, BH or "all"
66  timeClname(ps.getParameter<std::string>("timeClname")),
67  timeOffset(ps.getParameter<double>("timeOffset")),
68  nHitsTime(ps.getParameter<unsigned int>("nHitsTime")) {
69  if (detector == "HFNose") {
70  hits_hfnose_token = consumes<HGCRecHitCollection>(ps.getParameter<edm::InputTag>("HFNoseInput"));
72  } else if (detector == "all") {
73  hits_ee_token = consumes<HGCRecHitCollection>(ps.getParameter<edm::InputTag>("HGCEEInput"));
74  hits_fh_token = consumes<HGCRecHitCollection>(ps.getParameter<edm::InputTag>("HGCFHInput"));
75  hits_bh_token = consumes<HGCRecHitCollection>(ps.getParameter<edm::InputTag>("HGCBHInput"));
77  } else if (detector == "EE") {
78  hits_ee_token = consumes<HGCRecHitCollection>(ps.getParameter<edm::InputTag>("HGCEEInput"));
80  } else if (detector == "FH") {
81  hits_fh_token = consumes<HGCRecHitCollection>(ps.getParameter<edm::InputTag>("HGCFHInput"));
83  } else {
84  hits_bh_token = consumes<HGCRecHitCollection>(ps.getParameter<edm::InputTag>("HGCBHInput"));
86  }
87 
88  auto pluginPSet = ps.getParameter<edm::ParameterSet>("plugin");
89  if (detector == "HFNose") {
90  algo = HGCalLayerClusterAlgoFactory::get()->create("HFNoseCLUE", pluginPSet, consumesCollector());
91  algo->setAlgoId(algoId, true);
92  } else {
94  pluginPSet.getParameter<std::string>("type"), pluginPSet, consumesCollector());
95  algo->setAlgoId(algoId);
96  }
97 
98  produces<std::vector<float>>("InitialLayerClustersMask");
99  produces<std::vector<reco::BasicCluster>>();
100  produces<std::vector<reco::BasicCluster>>("sharing");
101  //density
102  produces<Density>();
103  //time for layer clusters
104  produces<edm::ValueMap<std::pair<float, float>>>(timeClname);
105 }
106 
108  // hgcalLayerClusters
110  edm::ParameterSetDescription pluginDesc;
111  pluginDesc.addNode(edm::PluginDescription<HGCalLayerClusterAlgoFactory>("type", "CLUE", true));
112 
113  desc.add<edm::ParameterSetDescription>("plugin", pluginDesc);
114  desc.add<std::string>("detector", "all")
115  ->setComment("all (does not include HFNose); other options: EE, FH, HFNose; other value defaults to BH");
116  desc.add<bool>("doSharing", false);
117  desc.add<edm::InputTag>("HFNoseInput", edm::InputTag("HGCalRecHit", "HGCHFNoseRecHits"));
118  desc.add<edm::InputTag>("HGCEEInput", edm::InputTag("HGCalRecHit", "HGCEERecHits"));
119  desc.add<edm::InputTag>("HGCFHInput", edm::InputTag("HGCalRecHit", "HGCHEFRecHits"));
120  desc.add<edm::InputTag>("HGCBHInput", edm::InputTag("HGCalRecHit", "HGCHEBRecHits"));
121  desc.add<std::string>("timeClname", "timeLayerCluster");
122  desc.add<double>("timeOffset", 0.0);
123  desc.add<unsigned int>("nHitsTime", 3);
124  descriptions.add("hgcalLayerClusters", desc);
125 }
126 
132 
133  std::unique_ptr<std::vector<reco::BasicCluster>> clusters(new std::vector<reco::BasicCluster>),
134  clusters_sharing(new std::vector<reco::BasicCluster>);
135  auto density = std::make_unique<Density>();
136 
137  algo->getEventSetup(es);
138 
139  //make a map detid-rechit
140  // NB for the moment just host EE and FH hits
141  // timing in digi for BH not implemented for now
142  std::unordered_map<uint32_t, const HGCRecHit*> hitmap;
143 
144  switch (algoId) {
146  evt.getByToken(hits_hfnose_token, hfnose_hits);
147  algo->populate(*hfnose_hits);
148  for (auto const& it : *hfnose_hits)
149  hitmap[it.detid().rawId()] = &(it);
150  break;
152  evt.getByToken(hits_ee_token, ee_hits);
153  algo->populate(*ee_hits);
154  for (auto const& it : *ee_hits)
155  hitmap[it.detid().rawId()] = &(it);
156  break;
158  evt.getByToken(hits_fh_token, fh_hits);
159  evt.getByToken(hits_bh_token, bh_hits);
160  if (fh_hits.isValid()) {
161  algo->populate(*fh_hits);
162  for (auto const& it : *fh_hits)
163  hitmap[it.detid().rawId()] = &(it);
164  } else if (bh_hits.isValid()) {
165  algo->populate(*bh_hits);
166  }
167  break;
169  evt.getByToken(hits_ee_token, ee_hits);
170  algo->populate(*ee_hits);
171  for (auto const& it : *ee_hits) {
172  hitmap[it.detid().rawId()] = &(it);
173  }
174  evt.getByToken(hits_fh_token, fh_hits);
175  algo->populate(*fh_hits);
176  for (auto const& it : *fh_hits) {
177  hitmap[it.detid().rawId()] = &(it);
178  }
179  evt.getByToken(hits_bh_token, bh_hits);
180  algo->populate(*bh_hits);
181  break;
182  default:
183  break;
184  }
185  algo->makeClusters();
186  *clusters = algo->getClusters(false);
187  if (doSharing)
188  *clusters_sharing = algo->getClusters(true);
189 
190  auto clusterHandle = evt.put(std::move(clusters));
191  auto clusterHandleSharing = evt.put(std::move(clusters_sharing), "sharing");
192 
193  //Keep the density
194  *density = algo->getDensity();
195  evt.put(std::move(density));
196 
197  edm::PtrVector<reco::BasicCluster> clusterPtrs, clusterPtrsSharing;
198 
199  std::vector<std::pair<float, float>> times;
200  times.reserve(clusterHandle->size());
201 
202  for (unsigned i = 0; i < clusterHandle->size(); ++i) {
203  edm::Ptr<reco::BasicCluster> ptr(clusterHandle, i);
204  clusterPtrs.push_back(ptr);
205 
206  std::pair<float, float> timeCl(-99., -1.);
207 
208  const reco::CaloCluster& sCl = (*clusterHandle)[i];
209  if (sCl.size() >= nHitsTime) {
210  std::vector<float> timeClhits;
211  std::vector<float> timeErrorClhits;
212 
213  for (auto const& hit : sCl.hitsAndFractions()) {
214  auto finder = hitmap.find(hit.first);
215  if (finder == hitmap.end())
216  continue;
217 
218  //time is computed wrt 0-25ns + offset and set to -1 if no time
219  const HGCRecHit* rechit = finder->second;
220  float rhTimeE = rechit->timeError();
221  //check on timeError to exclude scintillator
222  if (rhTimeE < 0.)
223  continue;
224  timeClhits.push_back(rechit->time() - timeOffset);
225  timeErrorClhits.push_back(1. / (rhTimeE * rhTimeE));
226  }
228  timeCl = timeEstimator.fixSizeHighestDensity(timeClhits, timeErrorClhits, nHitsTime);
229  }
230  times.push_back(timeCl);
231  }
232  std::unique_ptr<std::vector<float>> layerClustersMask(new std::vector<float>);
233  layerClustersMask->resize(clusterHandle->size(), 1.0);
234  evt.put(std::move(layerClustersMask), "InitialLayerClustersMask");
235 
236  auto timeCl = std::make_unique<edm::ValueMap<std::pair<float, float>>>();
237  edm::ValueMap<std::pair<float, float>>::Filler filler(*timeCl);
238  filler.insert(clusterHandle, times.begin(), times.end());
239  filler.fill();
240  evt.put(std::move(timeCl), timeClname);
241 
242  if (doSharing) {
243  for (unsigned i = 0; i < clusterHandleSharing->size(); ++i) {
244  edm::Ptr<reco::BasicCluster> ptr(clusterHandleSharing, i);
245  clusterPtrsSharing.push_back(ptr);
246  }
247  }
248  algo->reset();
249 }
250 
251 #endif //__RecoLocalCalo_HGCRecProducers_HGCalLayerClusterProducer_H__
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
void produce(edm::Event &, const edm::EventSetup &) override
reco::CaloCluster::AlgoId algoId
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void push_back(Ptr< T > const &iPtr)
Definition: PtrVector.h:149
ParameterDescriptionNode * addNode(ParameterDescriptionNode const &node)
edm::EDGetTokenT< HGCRecHitCollection > hits_hfnose_token
std::unique_ptr< HGCalClusteringAlgoBase > algo
std::map< DetId, float > Density
edm::EDGetTokenT< HGCRecHitCollection > hits_bh_token
std::pair< float, float > fixSizeHighestDensity(std::vector< float > &time, std::vector< float > weight=std::vector< float >(), unsigned int minNhits=3, float deltaT=0.210, float timeWidthBy=0.5)
float timeError() const
Definition: HGCRecHit.cc:79
const std::vector< std::pair< DetId, float > > & hitsAndFractions() const
Definition: CaloCluster.h:210
def move
Definition: eostools.py:511
constexpr float time() const
Definition: CaloRecHit.h:31
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
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
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:95
size_t size() const
size in number of hits (e.g. in crystals for ECAL)
Definition: CaloCluster.h:187
edm::EDGetTokenT< HGCRecHitCollection > hits_fh_token
#define get
hgcal_clustering::Density Density