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_
 
const bool calculatePositionInAlgo_
 
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 41 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 115 of file HGCalLayerClusterProducer.cc.

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

117  detector_(ps.getParameter<std::string>("detector")), // one of EE, FH, BH, HFNose
118  timeClname_(ps.getParameter<std::string>("timeClname")),
119  hitsTime_(ps.getParameter<unsigned int>("nHitsTime")),
120  caloGeomToken_(consumesCollector().esConsumes<CaloGeometry, CaloGeometryRecord>()),
121  calculatePositionInAlgo_(ps.getParameter<bool>("calculatePositionInAlgo")) {
122 #if DEBUG_CLUSTERS_ALPAKA
123  moduleLabel_ = ps.getParameter<std::string>("@module_label");
124 #endif
125  setAlgoId(); //sets algo id according to detector type
126  hits_token_ = consumes<HGCRecHitCollection>(ps.getParameter<edm::InputTag>("recHits"));
127 
128  auto pluginPSet = ps.getParameter<edm::ParameterSet>("plugin");
129  if (detector_ == "HFNose") {
130  algo_ = HGCalLayerClusterAlgoFactory::get()->create("HFNoseCLUE", pluginPSet);
131  algo_->setAlgoId(algoId_, true);
132  } else {
133  algo_ = HGCalLayerClusterAlgoFactory::get()->create(pluginPSet.getParameter<std::string>("type"), pluginPSet);
134  algo_->setAlgoId(algoId_);
135  }
136  thresholdW0_ = pluginPSet.getParameter<std::vector<double>>("thresholdW0");
137  positionDeltaRho2_ = pluginPSet.getParameter<double>("positionDeltaRho2");
138 
139  produces<std::vector<float>>("InitialLayerClustersMask");
140  produces<std::vector<reco::BasicCluster>>();
141  //time for layer clusters
142  produces<edm::ValueMap<std::pair<float, float>>>(timeClname_);
143 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
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 51 of file HGCalLayerClusterProducer.cc.

51 {}

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 160 of file HGCalLayerClusterProducer.cc.

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

Referenced by produce().

162  {
163  float total_weight = 0.f;
164  float maxEnergyValue = 0.f;
165  DetId maxEnergyIndex;
166  float x = 0.f;
167  float y = 0.f;
168 
169  for (auto const& hit : hitsAndFractions) {
170  //time is computed wrt 0-25ns + offset and set to -1 if no time
171  const HGCRecHit* rechit = hitmap[hit.first];
172  total_weight += rechit->energy();
173  if (rechit->energy() > maxEnergyValue) {
174  maxEnergyValue = rechit->energy();
175  maxEnergyIndex = rechit->detid();
176  }
177  }
178  float total_weight_log = 0.f;
179  auto thick = rhtools_.getSiThickIndex(maxEnergyIndex);
180  const GlobalPoint positionMaxEnergy(rhtools_.getPosition(maxEnergyIndex));
181  for (auto const& hit : hitsAndFractions) {
182  //time is computed wrt 0-25ns + offset and set to -1 if no time
183  const HGCRecHit* rechit = hitmap[hit.first];
184 
185  const GlobalPoint position(rhtools_.getPosition(rechit->detid()));
186 
187  if (thick != -1) { //silicon
188  //for silicon only just use 1+6 cells = 1.3cm for all thicknesses
189  const float d1 = position.x() - positionMaxEnergy.x();
190  const float d2 = position.y() - positionMaxEnergy.y();
191  if ((d1 * d1 + d2 * d2) > positionDeltaRho2_)
192  continue;
193 
194  float Wi = std::max(thresholdW0_[thick] + std::log(rechit->energy() / total_weight), 0.);
195  x += position.x() * Wi;
196  y += position.y() * Wi;
197  total_weight_log += Wi;
198  } else { //scintillator
199  x += position.x() * rechit->energy();
200  y += position.y() * rechit->energy();
201  }
202  }
203  if (thick != -1) {
204  total_weight = total_weight_log;
205  }
206  if (total_weight != 0.) {
207  float inv_tot_weight = 1.f / total_weight;
208  return math::XYZPoint(x * inv_tot_weight, y * inv_tot_weight, positionMaxEnergy.z());
209  } else {
210  return {positionMaxEnergy.x(), positionMaxEnergy.y(), positionMaxEnergy.z()};
211  }
212 }
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:140
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:216

◆ 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 214 of file HGCalLayerClusterProducer.cc.

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

Referenced by produce().

217  {
218  std::pair<float, float> timeCl(-99., -1.);
219 
220  if (sizeCluster >= hitsTime_) {
221  std::vector<float> timeClhits;
222  std::vector<float> timeErrorClhits;
223 
224  for (auto const& hit : hitsAndFractions) {
225  //time is computed wrt 0-25ns + offset and set to -1 if no time
226  const HGCRecHit* rechit = hitmap[hit.first];
227 
228  float rhTimeE = rechit->timeError();
229  //check on timeError to exclude scintillator
230  if (rhTimeE < 0.f)
231  continue;
232  timeClhits.push_back(rechit->time());
233  timeErrorClhits.push_back(1.f / (rhTimeE * rhTimeE));
234  }
236  timeCl = timeEstimator.fixSizeHighestDensity(timeClhits, timeErrorClhits, hitsTime_);
237  }
238  return timeCl;
239 }
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)
double f[11][100]
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 145 of file HGCalLayerClusterProducer.cc.

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

145  {
146  // hgcalLayerClusters
148  edm::ParameterSetDescription pluginDesc;
149  pluginDesc.addNode(edm::PluginDescription<HGCalLayerClusterAlgoFactory>("type", "SiCLUE", true));
150 
151  desc.add<edm::ParameterSetDescription>("plugin", pluginDesc);
152  desc.add<std::string>("detector", "EE")->setComment("options EE, FH, BH, HFNose; other value defaults to EE");
153  desc.add<edm::InputTag>("recHits", edm::InputTag("HGCalRecHit", "HGCEERecHits"));
154  desc.add<std::string>("timeClname", "timeLayerCluster");
155  desc.add<unsigned int>("nHitsTime", 3);
156  desc.add<bool>("calculatePositionInAlgo", true);
157  descriptions.add("hgcalLayerClusters", desc);
158 }
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 240 of file HGCalLayerClusterProducer.cc.

References algo_, calculatePosition(), calculatePositionInAlgo_, calculateTime(), caloGeomToken_, bsc_activity_cfg::clusters, detector_, hgcalUtils::DumpClusters::dumpInfos(), edm::EventID::event(), edm::Event::eventAuxiliary(), f, trigObjTnPSource_cfi::filler, relativeConstraints::geom, edm::Event::getByToken(), edm::EventSetup::getHandle(), hfClusterShapes_cfi::hits, hits_token_, reco::CaloCluster::hitsAndFractions(), mps_fire::i, edm::EventAuxiliary::id(), ALPAKA_ACCELERATOR_NAMESPACE::vertexFinder::it, edm::EventAuxiliary::luminosityBlock(), eostools::move(), edm::Event::put(), rhtools_, edm::EventAuxiliary::run(), convertSQLiteXML::runNumber, hgcal::RecHitTools::setGeometry(), reco::CaloCluster::setPosition(), reco::CaloCluster::size(), and timeClname_.

240  {
242 
243  std::unique_ptr<std::vector<reco::BasicCluster>> clusters(new std::vector<reco::BasicCluster>);
244 
247  algo_->getEventSetup(es, rhtools_);
248 
249  //make a map detid-rechit
250  // NB for the moment just host EE and FH hits
251  // timing in digi for BH not implemented for now
252  std::unordered_map<uint32_t, const HGCRecHit*> hitmap;
253 
255  algo_->populate(*hits);
256  for (auto const& it : *hits) {
257  hitmap[it.detid().rawId()] = &(it);
258  }
259 
260  algo_->makeClusters();
261  *clusters = algo_->getClusters(false);
262 
263  std::vector<std::pair<float, float>> times;
264  times.reserve(clusters->size());
265 
266  for (unsigned i = 0; i < clusters->size(); ++i) {
267  reco::CaloCluster& sCl = (*clusters)[i];
269  sCl.setPosition(calculatePosition(hitmap, sCl.hitsAndFractions()));
270  }
271  if (detector_ != "BH") {
272  times.push_back(calculateTime(hitmap, sCl.hitsAndFractions(), sCl.size()));
273  } else {
274  times.push_back(std::pair<float, float>(-99.f, -1.f));
275  }
276  }
277 
278 #if DEBUG_CLUSTERS_ALPAKA
280  auto runNumber = evt.eventAuxiliary().run();
281  auto lumiNumber = evt.eventAuxiliary().luminosityBlock();
282  auto evtNumber = evt.eventAuxiliary().id().event();
283 
284  dumper.dumpInfos(*clusters, moduleLabel_, runNumber, lumiNumber, evtNumber, true);
285 #endif
286 
287  auto clusterHandle = evt.put(std::move(clusters));
288 
289  if (detector_ == "HFNose") {
290  std::unique_ptr<std::vector<float>> layerClustersMask(new std::vector<float>);
291  layerClustersMask->resize(clusterHandle->size(), 1.0);
292  evt.put(std::move(layerClustersMask), "InitialLayerClustersMask");
293  }
294 
295  auto timeCl = std::make_unique<edm::ValueMap<std::pair<float, float>>>();
297  filler.insert(clusterHandle, times.begin(), times.end());
298  filler.fill();
299  evt.put(std::move(timeCl), timeClname_);
300 
301  algo_->reset();
302 }
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:210
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.
LuminosityBlockNumber_t luminosityBlock() const
void setPosition(const math::XYZPoint &p)
Definition: CaloCluster.h:140
edm::ESGetToken< CaloGeometry, CaloGeometryRecord > caloGeomToken_
void dumpInfos(const T &clusters, const std::string &moduleLabel, edm::RunNumber_t run, edm::LuminosityBlockNumber_t lumi, edm::EventNumber_t event, bool dumpCellsDetId=false) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:528
EventAuxiliary const & eventAuxiliary() const override
Definition: Event.h:95
EventID const & id() const
std::unique_ptr< HGCalClusteringAlgoBase > algo_
double f[11][100]
size_t size() const
size in number of hits (e.g. in crystals for ECAL)
Definition: CaloCluster.h:187
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:79
def move(src, dest)
Definition: eostools.py:511
EventNumber_t event() const
Definition: EventID.h:40
RunNumber_t run() const

◆ setAlgoId()

void HGCalLayerClusterProducer::setAlgoId ( )
private

Member Data Documentation

◆ algo_

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

Definition at line 72 of file HGCalLayerClusterProducer.cc.

Referenced by HGCalLayerClusterProducer(), and produce().

◆ algoId_

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

Definition at line 70 of file HGCalLayerClusterProducer.cc.

Referenced by HGCalLayerClusterProducer(), and setAlgoId().

◆ calculatePositionInAlgo_

const bool HGCalLayerClusterProducer::calculatePositionInAlgo_
private

Definition at line 83 of file HGCalLayerClusterProducer.cc.

Referenced by produce().

◆ caloGeomToken_

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

Definition at line 82 of file HGCalLayerClusterProducer.cc.

Referenced by produce().

◆ detector_

std::string HGCalLayerClusterProducer::detector_
private

Definition at line 73 of file HGCalLayerClusterProducer.cc.

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

◆ hits_token_

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

Definition at line 68 of file HGCalLayerClusterProducer.cc.

Referenced by HGCalLayerClusterProducer(), and produce().

◆ hitsTime_

unsigned int HGCalLayerClusterProducer::hitsTime_
private

Definition at line 76 of file HGCalLayerClusterProducer.cc.

Referenced by calculateTime().

◆ positionDeltaRho2_

double HGCalLayerClusterProducer::positionDeltaRho2_
private

Definition at line 80 of file HGCalLayerClusterProducer.cc.

Referenced by calculatePosition(), and HGCalLayerClusterProducer().

◆ rhtools_

hgcal::RecHitTools HGCalLayerClusterProducer::rhtools_
private

Definition at line 81 of file HGCalLayerClusterProducer.cc.

Referenced by calculatePosition(), and produce().

◆ thresholdW0_

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

Definition at line 79 of file HGCalLayerClusterProducer.cc.

Referenced by calculatePosition(), and HGCalLayerClusterProducer().

◆ timeClname_

std::string HGCalLayerClusterProducer::timeClname_
private

Definition at line 75 of file HGCalLayerClusterProducer.cc.

Referenced by HGCalLayerClusterProducer(), and produce().