CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Private Attributes | Static Private Attributes
HGCalDepthPreClusterer Class Reference

#include <HGCalDepthPreClusterer.h>

Public Types

typedef std::vector< reco::BasicClusterClusterCollection
 

Public Member Functions

void getEvent (const edm::Event &ev)
 
void getEventSetup (const edm::EventSetup &es)
 
 HGCalDepthPreClusterer ()
 
 HGCalDepthPreClusterer (const edm::ParameterSet &conf, edm::ConsumesCollector &sumes, std::vector< float > radii_in, uint32_t min_clusters, bool real_space_cone)
 
std::vector< reco::HGCalMultiClustermakePreClusters (const reco::HGCalMultiCluster::ClusterCollection &) const
 

Private Attributes

std::unique_ptr< hgcal::ClusterToolsclusterTools
 
uint32_t minClusters
 
std::vector< float > radii
 
bool realSpaceCone
 

Static Private Attributes

static const unsigned int lastLayerBH = 52
 
static const unsigned int lastLayerEE = 28
 
static const unsigned int lastLayerFH = 40
 

Detailed Description

Definition at line 13 of file HGCalDepthPreClusterer.h.

Member Typedef Documentation

Definition at line 31 of file HGCalDepthPreClusterer.h.

Constructor & Destructor Documentation

HGCalDepthPreClusterer::HGCalDepthPreClusterer ( )
inline

Definition at line 17 of file HGCalDepthPreClusterer.h.

References clusterTools, minClusters, and realSpaceCone.

17  : radii({0.,0.,0.,}), minClusters(0), realSpaceCone(false), clusterTools(nullptr)
18  {
19  }
std::unique_ptr< hgcal::ClusterTools > clusterTools
std::vector< float > radii
HGCalDepthPreClusterer::HGCalDepthPreClusterer ( const edm::ParameterSet conf,
edm::ConsumesCollector sumes,
std::vector< float >  radii_in,
uint32_t  min_clusters,
bool  real_space_cone 
)
inline

Definition at line 21 of file HGCalDepthPreClusterer.h.

21  :
22  radii(radii_in),
23  minClusters(min_clusters),
24  realSpaceCone(real_space_cone),
25  clusterTools(std::make_unique<hgcal::ClusterTools>(conf,sumes)) {
26  }
std::unique_ptr< hgcal::ClusterTools > clusterTools
std::vector< float > radii

Member Function Documentation

void HGCalDepthPreClusterer::getEvent ( const edm::Event ev)
inline

Definition at line 28 of file HGCalDepthPreClusterer.h.

References clusterTools.

28 { clusterTools->getEvent(ev); }
std::unique_ptr< hgcal::ClusterTools > clusterTools
void HGCalDepthPreClusterer::getEventSetup ( const edm::EventSetup es)
inline

Definition at line 29 of file HGCalDepthPreClusterer.h.

References clusterTools.

29 { clusterTools->getEventSetup(es); }
std::unique_ptr< hgcal::ClusterTools > clusterTools
std::vector< reco::HGCalMultiCluster > HGCalDepthPreClusterer::makePreClusters ( const reco::HGCalMultiCluster::ClusterCollection thecls) const

Definition at line 36 of file HGCalDepthPreClusterer.cc.

References clusterTools, mps_fire::i, lastLayerBH, lastLayerEE, lastLayerFH, minClusters, reco::HGCalMultiCluster::push_back(), radii, realSpaceCone, reco::HGCalMultiCluster::size(), sorted_indices(), groupFilesInBlocks::temp, and z.

36  {
37 
38  std::vector<reco::HGCalMultiCluster> thePreClusters;
39  std::vector<size_t> es = sorted_indices(thecls);
40  std::vector<int> vused(es.size(),0);
41  unsigned int used = 0;
42 
43  for(unsigned int i = 0; i < es.size(); ++i) {
44  if(vused[i]==0) {
46  temp.push_back(thecls[es[i]]);
47  vused[i]=(thecls[es[i]]->z()>0)? 1 : -1;
48  ++used;
49  for(unsigned int j = i+1; j < es.size(); ++j) {
50  if(vused[j]==0) {
51  float distanceCheck = 9999.;
52  if( realSpaceCone ) distanceCheck = distAxisCluster2(thecls[es[i]],thecls[es[j]]);
53  else distanceCheck = dist2(thecls[es[i]],thecls[es[j]]);
54  DetId detid = thecls[es[j]]->hitsAndFractions()[0].first();
55  unsigned int layer = clusterTools->getLayer(detid);
56  float radius2 = 9999.;
57  if(layer <= lastLayerEE) radius2 = radii[0]*radii[0];
58  else if(layer <= lastLayerFH) radius2 = radii[1]*radii[1];
59  else if(layer <= lastLayerBH) radius2 = radii[2]*radii[2];
60  else assert(radius2<100. && "nonsense layer value - cannot assign multicluster radius");
61  if( distanceCheck<radius2 && int(thecls[es[j]]->z()*vused[i])>0 ) {
62  temp.push_back(thecls[es[j]]);
63  vused[j]=vused[i];
64  ++used;
65  }
66  }
67  }
68  if( temp.size() > minClusters ) {
69  thePreClusters.push_back(temp);
70  auto& back = thePreClusters.back();
71  back.setPosition(clusterTools->getMultiClusterPosition(back));
72  back.setEnergy(clusterTools->getMultiClusterEnergy(back));
73  }
74  }
75  }
76 
77 
78 
79  return thePreClusters;
80 }
static const unsigned int lastLayerBH
std::unique_ptr< hgcal::ClusterTools > clusterTools
unsigned int size() const
void push_back(const edm::Ptr< reco::BasicCluster > &b)
std::vector< size_t > sorted_indices(const std::vector< T > &v)
Definition: DetId.h:18
std::vector< float > radii
static const unsigned int lastLayerEE
static const unsigned int lastLayerFH

Member Data Documentation

std::unique_ptr<hgcal::ClusterTools> HGCalDepthPreClusterer::clusterTools
private
const unsigned int HGCalDepthPreClusterer::lastLayerBH = 52
staticprivate

Definition at line 45 of file HGCalDepthPreClusterer.h.

Referenced by makePreClusters().

const unsigned int HGCalDepthPreClusterer::lastLayerEE = 28
staticprivate

Definition at line 43 of file HGCalDepthPreClusterer.h.

Referenced by makePreClusters().

const unsigned int HGCalDepthPreClusterer::lastLayerFH = 40
staticprivate

Definition at line 44 of file HGCalDepthPreClusterer.h.

Referenced by makePreClusters().

uint32_t HGCalDepthPreClusterer::minClusters
private

Definition at line 38 of file HGCalDepthPreClusterer.h.

Referenced by HGCalDepthPreClusterer(), and makePreClusters().

std::vector<float> HGCalDepthPreClusterer::radii
private

Definition at line 37 of file HGCalDepthPreClusterer.h.

Referenced by makePreClusters().

bool HGCalDepthPreClusterer::realSpaceCone
private

flag to use cartesian space clustering.

Definition at line 39 of file HGCalDepthPreClusterer.h.

Referenced by HGCalDepthPreClusterer(), and makePreClusters().