CMS 3D CMS Logo

HGCalLayerClusterProducer.cc
Go to the documentation of this file.
1 // Authors: Olivie Franklova - olivie.abigail.franklova@cern.ch
2 // Date: 03/2023
3 // @file create layer clusters
4 
5 // user include files
9 
15 
23 
26 
29 
32 
34 
36 public:
51  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
52 
59  void produce(edm::Event&, const edm::EventSetup&) override;
60 
61 private:
63 
65 
66  std::unique_ptr<HGCalClusteringAlgoBase> algo_;
68 
70  unsigned int hitsTime_;
71 
72  // for calculate position
73  std::vector<double> thresholdW0_;
77 
81  void setAlgoId();
82 
90  math::XYZPoint calculatePosition(std::unordered_map<uint32_t, const HGCRecHit*>& hitmap,
91  const std::vector<std::pair<DetId, float>>& hitsAndFractions);
92 
100  std::pair<float, float> calculateTime(std::unordered_map<uint32_t, const HGCRecHit*>& hitmap,
101  const std::vector<std::pair<DetId, float>>& hitsAndFractions,
102  size_t sizeCluster);
103 };
104 
106  : algoId_(reco::CaloCluster::undefined),
107  detector_(ps.getParameter<std::string>("detector")), // one of EE, FH, BH, HFNose
108  timeClname_(ps.getParameter<std::string>("timeClname")),
109  hitsTime_(ps.getParameter<unsigned int>("nHitsTime")),
110  caloGeomToken_(consumesCollector().esConsumes<CaloGeometry, CaloGeometryRecord>()) {
111  setAlgoId(); //sets algo id according to detector type
112  hits_token_ = consumes<HGCRecHitCollection>(ps.getParameter<edm::InputTag>("recHits"));
113 
114  auto pluginPSet = ps.getParameter<edm::ParameterSet>("plugin");
115  if (detector_ == "HFNose") {
116  algo_ = HGCalLayerClusterAlgoFactory::get()->create("HFNoseCLUE", pluginPSet);
117  algo_->setAlgoId(algoId_, true);
118  } else {
119  algo_ = HGCalLayerClusterAlgoFactory::get()->create(pluginPSet.getParameter<std::string>("type"), pluginPSet);
120  algo_->setAlgoId(algoId_);
121  }
122  thresholdW0_ = pluginPSet.getParameter<std::vector<double>>("thresholdW0");
123  positionDeltaRho2_ = pluginPSet.getParameter<double>("positionDeltaRho2");
124 
125  produces<std::vector<float>>("InitialLayerClustersMask");
126  produces<std::vector<reco::BasicCluster>>();
127  //time for layer clusters
128  produces<edm::ValueMap<std::pair<float, float>>>(timeClname_);
129 }
130 
132  // hgcalLayerClusters
134  edm::ParameterSetDescription pluginDesc;
135  pluginDesc.addNode(edm::PluginDescription<HGCalLayerClusterAlgoFactory>("type", "SiCLUE", true));
136 
137  desc.add<edm::ParameterSetDescription>("plugin", pluginDesc);
138  desc.add<std::string>("detector", "EE")->setComment("options EE, FH, BH, HFNose; other value defaults to EE");
139  desc.add<edm::InputTag>("recHits", edm::InputTag("HGCalRecHit", "HGCEERecHits"));
140  desc.add<std::string>("timeClname", "timeLayerCluster");
141  desc.add<unsigned int>("nHitsTime", 3);
142  descriptions.add("hgcalLayerClusters", desc);
143 }
144 
146  std::unordered_map<uint32_t, const HGCRecHit*>& hitmap,
147  const std::vector<std::pair<DetId, float>>& hitsAndFractions) {
148  float total_weight = 0.f;
149  float maxEnergyValue = 0.f;
150  DetId maxEnergyIndex;
151  float x = 0.f;
152  float y = 0.f;
153 
154  for (auto const& hit : hitsAndFractions) {
155  //time is computed wrt 0-25ns + offset and set to -1 if no time
156  const HGCRecHit* rechit = hitmap[hit.first];
157  total_weight += rechit->energy();
158  if (rechit->energy() > maxEnergyValue) {
159  maxEnergyValue = rechit->energy();
160  maxEnergyIndex = rechit->detid();
161  }
162  }
163  float total_weight_log = 0.f;
164  auto thick = rhtools_.getSiThickIndex(maxEnergyIndex);
165  const GlobalPoint positionMaxEnergy(rhtools_.getPosition(maxEnergyIndex));
166  for (auto const& hit : hitsAndFractions) {
167  //time is computed wrt 0-25ns + offset and set to -1 if no time
168  const HGCRecHit* rechit = hitmap[hit.first];
169 
170  const GlobalPoint position(rhtools_.getPosition(rechit->detid()));
171 
172  if (thick != -1) { //silicon
173  //for silicon only just use 1+6 cells = 1.3cm for all thicknesses
174  const float d1 = position.x() - positionMaxEnergy.x();
175  const float d2 = position.y() - positionMaxEnergy.y();
176  if ((d1 * d1 + d2 * d2) > positionDeltaRho2_)
177  continue;
178 
179  float Wi = std::max(thresholdW0_[thick] + std::log(rechit->energy() / total_weight), 0.);
180  x += position.x() * Wi;
181  y += position.y() * Wi;
182  total_weight_log += Wi;
183  } else { //scintillator
184  x += position.x() * rechit->energy();
185  y += position.y() * rechit->energy();
186  }
187  }
188  if (thick != -1) {
189  total_weight = total_weight_log;
190  }
191  if (total_weight != 0.) {
192  float inv_tot_weight = 1.f / total_weight;
193  return math::XYZPoint(x * inv_tot_weight, y * inv_tot_weight, positionMaxEnergy.z());
194  } else {
195  return math::XYZPoint(0.f, 0.f, 0.f);
196  }
197 }
198 
200  std::unordered_map<uint32_t, const HGCRecHit*>& hitmap,
201  const std::vector<std::pair<DetId, float>>& hitsAndFractions,
202  size_t sizeCluster) {
203  std::pair<float, float> timeCl(-99., -1.);
204 
205  if (sizeCluster >= hitsTime_) {
206  std::vector<float> timeClhits;
207  std::vector<float> timeErrorClhits;
208 
209  for (auto const& hit : hitsAndFractions) {
210  //time is computed wrt 0-25ns + offset and set to -1 if no time
211  const HGCRecHit* rechit = hitmap[hit.first];
212 
213  float rhTimeE = rechit->timeError();
214  //check on timeError to exclude scintillator
215  if (rhTimeE < 0.)
216  continue;
217  timeClhits.push_back(rechit->time());
218  timeErrorClhits.push_back(1. / (rhTimeE * rhTimeE));
219  }
221  timeCl = timeEstimator.fixSizeHighestDensity(timeClhits, timeErrorClhits, hitsTime_);
222  }
223  return timeCl;
224 }
227 
228  std::unique_ptr<std::vector<reco::BasicCluster>> clusters(new std::vector<reco::BasicCluster>);
229 
232  algo_->getEventSetup(es, rhtools_);
233 
234  //make a map detid-rechit
235  // NB for the moment just host EE and FH hits
236  // timing in digi for BH not implemented for now
237  std::unordered_map<uint32_t, const HGCRecHit*> hitmap;
238 
240  algo_->populate(*hits);
241  for (auto const& it : *hits) {
242  hitmap[it.detid().rawId()] = &(it);
243  }
244 
245  algo_->makeClusters();
246  *clusters = algo_->getClusters(false);
247 
248  std::vector<std::pair<float, float>> times;
249  times.reserve(clusters->size());
250 
251  for (unsigned i = 0; i < clusters->size(); ++i) {
252  const reco::CaloCluster& sCl = (*clusters)[i];
253  (*clusters)[i].setPosition(calculatePosition(hitmap, sCl.hitsAndFractions()));
254  if (detector_ != "BH") {
255  times.push_back(calculateTime(hitmap, sCl.hitsAndFractions(), sCl.size()));
256  } else {
257  times.push_back(std::pair<float, float>(-99., -1.));
258  }
259  }
260 
261  auto clusterHandle = evt.put(std::move(clusters));
262 
263  if (detector_ == "HFNose") {
264  std::unique_ptr<std::vector<float>> layerClustersMask(new std::vector<float>);
265  layerClustersMask->resize(clusterHandle->size(), 1.0);
266  evt.put(std::move(layerClustersMask), "InitialLayerClustersMask");
267  }
268 
269  auto timeCl = std::make_unique<edm::ValueMap<std::pair<float, float>>>();
271  filler.insert(clusterHandle, times.begin(), times.end());
272  filler.fill();
273  evt.put(std::move(timeCl), timeClname_);
274 
275  algo_->reset();
276 }
277 
279  if (detector_ == "HFNose") {
281  } else if (detector_ == "EE") {
283  } else { //for FH or BH
285  }
286 }
287 
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
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
math::XYZPoint calculatePosition(std::unordered_map< uint32_t, const HGCRecHit *> &hitmap, const std::vector< std::pair< DetId, float >> &hitsAndFractions)
Counts position for all points in the cluster.
const std::vector< std::pair< DetId, float > > & hitsAndFractions() const
Definition: CaloCluster.h:209
void produce(edm::Event &, const edm::EventSetup &) override
Method run the algoritm to get clusters.
constexpr const DetId & detid() const
Definition: CaloRecHit.h:33
std::pair< float, float > calculateTime(std::unordered_map< uint32_t, const HGCRecHit *> &hitmap, const std::vector< std::pair< DetId, float >> &hitsAndFractions, size_t sizeCluster)
Counts time for all points in the cluster.
void setPosition(const math::XYZPoint &p)
Definition: CaloCluster.h:139
ParameterDescriptionNode * addNode(ParameterDescriptionNode const &node)
edm::ESGetToken< CaloGeometry, CaloGeometryRecord > caloGeomToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:528
constexpr float energy() const
Definition: CaloRecHit.h:29
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)
std::unique_ptr< HGCalClusteringAlgoBase > algo_
GlobalPoint getPosition(const DetId &id) const
Definition: RecHitTools.cc:140
reco::CaloCluster::AlgoId algoId_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Method fill description which will be used in pyhton file.
double f[11][100]
size_t size() const
size in number of hits (e.g. in crystals for ECAL)
Definition: CaloCluster.h:186
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void setAlgoId()
Sets algoId accordingly to the detector type.
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
Definition: DetId.h:17
edm::EDGetTokenT< HGCRecHitCollection > hits_token_
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
void add(std::string const &label, ParameterSetDescription const &psetDescription)
void setGeometry(CaloGeometry const &)
Definition: RecHitTools.cc:79
constexpr float time() const
Definition: CaloRecHit.h:31
HGCalLayerClusterProducer(const edm::ParameterSet &)
Constructor with parameter settings - which can be changed in hgcalLayerCluster_cff.py. Constructor will set all variables by input param ps. algoID variables will be set accordingly to the detector type.
fixed size matrix
static int position[264][3]
Definition: ReadPGInfo.cc:289
float timeError() const
Definition: HGCRecHit.cc:79
#define get
static constexpr float d1
int getSiThickIndex(const DetId &) const
Definition: RecHitTools.cc:216
def move(src, dest)
Definition: eostools.py:511