CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
HGCalLayerClusterProducer Class Reference
Inheritance diagram for HGCalLayerClusterProducer:
edm::stream::EDProducer<>

Public Member Functions

 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. More...
 
void produce (edm::Event &, const edm::EventSetup &) override
 Method run the algoritm to get clusters. More...
 
 ~HGCalLayerClusterProducer () override
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 Method fill description which will be used in pyhton file. More...
 

Private Member Functions

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. More...
 
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. More...
 
void setAlgoId ()
 Sets algoId accordingly to the detector type. More...
 

Private Attributes

std::unique_ptr< HGCalClusteringAlgoBasealgo_
 
reco::CaloCluster::AlgoId algoId_
 
edm::ESGetToken< CaloGeometry, CaloGeometryRecordcaloGeomToken_
 
std::string detector_
 
edm::EDGetTokenT< HGCRecHitCollectionhits_token_
 
unsigned int hitsTime_
 
double positionDeltaRho2_
 
hgcal::RecHitTools rhtools_
 
std::vector< double > thresholdW0_
 
std::string timeClname_
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 

Detailed Description

Definition at line 35 of file HGCalLayerClusterProducer.cc.

Constructor & Destructor Documentation

◆ HGCalLayerClusterProducer()

HGCalLayerClusterProducer::HGCalLayerClusterProducer ( const edm::ParameterSet ps)

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.

Parameters
[in]psparametr set to set variables

Definition at line 105 of file HGCalLayerClusterProducer.cc.

References algo_, algoId_, detector_, get, edm::ParameterSet::getParameter(), hits_token_, positionDeltaRho2_, setAlgoId(), AlCaHLTBitMon_QueryRunRegistry::string, thresholdW0_, and timeClname_.

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 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::ESGetToken< CaloGeometry, CaloGeometryRecord > caloGeomToken_
std::unique_ptr< HGCalClusteringAlgoBase > algo_
reco::CaloCluster::AlgoId algoId_
void setAlgoId()
Sets algoId accordingly to the detector type.
edm::EDGetTokenT< HGCRecHitCollection > hits_token_
#define get

◆ ~HGCalLayerClusterProducer()

HGCalLayerClusterProducer::~HGCalLayerClusterProducer ( )
inlineoverride

Definition at line 45 of file HGCalLayerClusterProducer.cc.

45 {}

Member Function Documentation

◆ calculatePosition()

math::XYZPoint HGCalLayerClusterProducer::calculatePosition ( std::unordered_map< uint32_t, const HGCRecHit *> &  hitmap,
const std::vector< std::pair< DetId, float >> &  hitsAndFractions 
)
private

Counts position for all points in the cluster.

Parameters
[in]hitmaphitmap to find correct RecHit
[in]hitsAndFractionall hits in the cluster
Returns
counted position

Definition at line 145 of file HGCalLayerClusterProducer.cc.

References d1, CaloRecHit::detid(), CaloRecHit::energy(), f, hgcal::RecHitTools::getPosition(), hgcal::RecHitTools::getSiThickIndex(), dqm-mbProfile::log, SiStripPI::max, position, positionDeltaRho2_, rhtools_, thresholdW0_, x, and y.

Referenced by produce().

147  {
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 }
constexpr const DetId & detid() const
Definition: CaloRecHit.h:33
constexpr float energy() const
Definition: CaloRecHit.h:29
GlobalPoint getPosition(const DetId &id) const
Definition: RecHitTools.cc:129
double f[11][100]
Definition: DetId.h:17
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
static int position[264][3]
Definition: ReadPGInfo.cc:289
static constexpr float d1
int getSiThickIndex(const DetId &) const
Definition: RecHitTools.cc:205

◆ calculateTime()

std::pair< float, float > HGCalLayerClusterProducer::calculateTime ( std::unordered_map< uint32_t, const HGCRecHit *> &  hitmap,
const std::vector< std::pair< DetId, float >> &  hitsAndFractions,
size_t  sizeCluster 
)
private

Counts time for all points in the cluster.

Parameters
[in]hitmaphitmap to find correct RecHit only for silicon (not for BH-HSci)
[in]hitsAndFractionall hits in the cluster
Returns
counted time

Definition at line 199 of file HGCalLayerClusterProducer.cc.

References hgcalsimclustertime::ComputeClusterTime::fixSizeHighestDensity(), hitsTime_, CaloRecHit::time(), and HGCRecHit::timeError().

Referenced by produce().

202  {
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 }
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)
constexpr float time() const
Definition: CaloRecHit.h:31
float timeError() const
Definition: HGCRecHit.cc:79

◆ fillDescriptions()

void HGCalLayerClusterProducer::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Method fill description which will be used in pyhton file.

Parameters
[out]descriptionto be fill

Definition at line 131 of file HGCalLayerClusterProducer.cc.

References edm::ConfigurationDescriptions::add(), edm::ParameterSetDescription::addNode(), submitPVResolutionJobs::desc, ProducerED_cfi::InputTag, and AlCaHLTBitMon_QueryRunRegistry::string.

131  {
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 }
ParameterDescriptionNode * addNode(ParameterDescriptionNode const &node)
void add(std::string const &label, ParameterSetDescription const &psetDescription)

◆ produce()

void HGCalLayerClusterProducer::produce ( edm::Event evt,
const edm::EventSetup es 
)
override

Method run the algoritm to get clusters.

Parameters
[in,out]evtfrom get info and put result
[in]esto get event setup info

Definition at line 225 of file HGCalLayerClusterProducer.cc.

References algo_, calculatePosition(), calculateTime(), caloGeomToken_, bsc_activity_cfg::clusters, detector_, trigObjTnPSource_cfi::filler, relativeConstraints::geom, edm::Event::getByToken(), edm::EventSetup::getHandle(), hfClusterShapes_cfi::hits, hits_token_, reco::CaloCluster::hitsAndFractions(), mps_fire::i, eostools::move(), edm::Event::put(), rhtools_, hgcal::RecHitTools::setGeometry(), reco::CaloCluster::setPosition(), reco::CaloCluster::size(), and timeClname_.

225  {
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 }
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
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
edm::ESGetToken< CaloGeometry, CaloGeometryRecord > caloGeomToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:536
std::unique_ptr< HGCalClusteringAlgoBase > algo_
size_t size() const
size in number of hits (e.g. in crystals for ECAL)
Definition: CaloCluster.h:186
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
edm::EDGetTokenT< HGCRecHitCollection > hits_token_
void setGeometry(CaloGeometry const &)
Definition: RecHitTools.cc:68
def move(src, dest)
Definition: eostools.py:511

◆ setAlgoId()

void HGCalLayerClusterProducer::setAlgoId ( )
private

Sets algoId accordingly to the detector type.

Definition at line 278 of file HGCalLayerClusterProducer.cc.

References algoId_, detector_, reco::CaloCluster::hfnose, reco::CaloCluster::hgcal_em, and reco::CaloCluster::hgcal_had.

Referenced by HGCalLayerClusterProducer().

278  {
279  if (detector_ == "HFNose") {
281  } else if (detector_ == "EE") {
283  } else { //for FH or BH
285  }
286 }
reco::CaloCluster::AlgoId algoId_

Member Data Documentation

◆ algo_

std::unique_ptr<HGCalClusteringAlgoBase> HGCalLayerClusterProducer::algo_
private

Definition at line 66 of file HGCalLayerClusterProducer.cc.

Referenced by HGCalLayerClusterProducer(), and produce().

◆ algoId_

reco::CaloCluster::AlgoId HGCalLayerClusterProducer::algoId_
private

Definition at line 64 of file HGCalLayerClusterProducer.cc.

Referenced by HGCalLayerClusterProducer(), and setAlgoId().

◆ caloGeomToken_

edm::ESGetToken<CaloGeometry, CaloGeometryRecord> HGCalLayerClusterProducer::caloGeomToken_
private

Definition at line 76 of file HGCalLayerClusterProducer.cc.

Referenced by produce().

◆ detector_

std::string HGCalLayerClusterProducer::detector_
private

Definition at line 67 of file HGCalLayerClusterProducer.cc.

Referenced by HGCalLayerClusterProducer(), produce(), and setAlgoId().

◆ hits_token_

edm::EDGetTokenT<HGCRecHitCollection> HGCalLayerClusterProducer::hits_token_
private

Definition at line 62 of file HGCalLayerClusterProducer.cc.

Referenced by HGCalLayerClusterProducer(), and produce().

◆ hitsTime_

unsigned int HGCalLayerClusterProducer::hitsTime_
private

Definition at line 70 of file HGCalLayerClusterProducer.cc.

Referenced by calculateTime().

◆ positionDeltaRho2_

double HGCalLayerClusterProducer::positionDeltaRho2_
private

Definition at line 74 of file HGCalLayerClusterProducer.cc.

Referenced by calculatePosition(), and HGCalLayerClusterProducer().

◆ rhtools_

hgcal::RecHitTools HGCalLayerClusterProducer::rhtools_
private

Definition at line 75 of file HGCalLayerClusterProducer.cc.

Referenced by calculatePosition(), and produce().

◆ thresholdW0_

std::vector<double> HGCalLayerClusterProducer::thresholdW0_
private

Definition at line 73 of file HGCalLayerClusterProducer.cc.

Referenced by calculatePosition(), and HGCalLayerClusterProducer().

◆ timeClname_

std::string HGCalLayerClusterProducer::timeClname_
private

Definition at line 69 of file HGCalLayerClusterProducer.cc.

Referenced by HGCalLayerClusterProducer(), and produce().