CMS 3D CMS Logo

List of all members | Classes | Public Types | Public Member Functions | Private Member Functions | Private Attributes
Multi5x5ClusterAlgo Class Reference

#include <Multi5x5ClusterAlgo.h>

Classes

class  ProtoBasicCluster
 

Public Types

typedef math::XYZPoint Point
 point in the space More...
 

Public Member Functions

std::vector< reco::BasicClustermakeClusters (const EcalRecHitCollection *hits, const CaloSubdetectorGeometry *geometry, const CaloSubdetectorTopology *topology_p, const CaloSubdetectorGeometry *geometryES_p, reco::CaloID::Detectors detector, bool regional=false, const std::vector< RectangularEtaPhiRegion > &regions=std::vector< RectangularEtaPhiRegion >())
 
 Multi5x5ClusterAlgo ()
 
 Multi5x5ClusterAlgo (double ebst, double ecst, const std::vector< int > &v_chstatus, const PositionCalc &posCalc, bool reassignSeedCrysToClusterItSeeds=false)
 
virtual ~Multi5x5ClusterAlgo ()
 

Private Member Functions

void addCrystal (const DetId &det)
 
bool checkMaxima (CaloNavigator< DetId > &navigator, const EcalRecHitCollection *hits)
 
void mainSearch (const EcalRecHitCollection *hits, const CaloSubdetectorGeometry *geometry_p, const CaloSubdetectorTopology *topology_p, const CaloSubdetectorGeometry *geometryES_p)
 
void makeCluster (const EcalRecHitCollection *hits, const CaloSubdetectorGeometry *geometry_p, const CaloSubdetectorGeometry *geometryES_p, const EcalRecHitCollection::const_iterator &seedIt, bool seedOutside)
 
void prepareCluster (CaloNavigator< DetId > &navigator, const EcalRecHitCollection *hits, const CaloSubdetectorGeometry *geometry)
 

Private Attributes

std::set< DetIdcanSeed_s
 
std::vector< reco::BasicClusterclusters_v
 
std::vector< std::pair< DetId, float > > current_v
 
reco::CaloID::Detectors detector_
 The ecal region used. More...
 
double ecalBarrelSeedThreshold
 
double ecalEndcapSeedThreshold
 
PositionCalc posCalculator_
 
std::vector< ProtoBasicClusterprotoClusters_
 
bool reassignSeedCrysToClusterItSeeds_
 
const EcalRecHitCollectionrecHits_
 
std::vector< EcalRecHitseeds
 
std::set< DetIdused_s
 
std::vector< int > v_chstatus_
 
std::vector< std::pair< DetId, int > > whichClusCrysBelongsTo_
 

Detailed Description

Definition at line 27 of file Multi5x5ClusterAlgo.h.

Member Typedef Documentation

point in the space

Definition at line 80 of file Multi5x5ClusterAlgo.h.

Constructor & Destructor Documentation

Multi5x5ClusterAlgo::Multi5x5ClusterAlgo ( )
inline

Definition at line 56 of file Multi5x5ClusterAlgo.h.

56  {
57  }
Multi5x5ClusterAlgo::Multi5x5ClusterAlgo ( double  ebst,
double  ecst,
const std::vector< int > &  v_chstatus,
const PositionCalc posCalc,
bool  reassignSeedCrysToClusterItSeeds = false 
)
inline
virtual Multi5x5ClusterAlgo::~Multi5x5ClusterAlgo ( )
inlinevirtual

Member Function Documentation

void Multi5x5ClusterAlgo::addCrystal ( const DetId det)
private

Definition at line 404 of file Multi5x5ClusterAlgo.cc.

405 {
406 
408  if ((thisIt != recHits_->end()) && (thisIt->id() != DetId(0)))
409  {
410  if ((used_s.find(thisIt->id()) == used_s.end()))
411  {
412  //std::cout << " ... this is a good crystal and will be added" << std::endl;
413  current_v.push_back( std::pair<DetId, float>(det, 1.) ); // by default hit energy fractions are set at 1.
414  used_s.insert(det);
415  }
416  }
417 
418 }
const EcalRecHitCollection * recHits_
std::vector< EcalRecHit >::const_iterator const_iterator
std::set< DetId > used_s
std::vector< std::pair< DetId, float > > current_v
const_iterator end() const
Definition: DetId.h:18
iterator find(key_type k)
bool Multi5x5ClusterAlgo::checkMaxima ( CaloNavigator< DetId > &  navigator,
const EcalRecHitCollection hits 
)
private

Definition at line 296 of file Multi5x5ClusterAlgo.cc.

References CaloNavigator< T, TOPO >::east(), spr::find(), edm::SortedCollection< T, SORT >::find(), CaloNavigator< T, TOPO >::home(), mps_fire::i, CaloNavigator< T, TOPO >::north(), CaloNavigator< T, TOPO >::pos(), CaloNavigator< T, TOPO >::south(), and CaloNavigator< T, TOPO >::west().

298 {
299 
300  bool maxima = true;
302  EcalRecHitCollection::const_iterator seedHit = hits->find(navigator.pos());
303  double seedEnergy = seedHit->energy();
304 
305  std::vector<DetId> swissCrossVec;
306  swissCrossVec.clear();
307 
308  swissCrossVec.push_back(navigator.west());
309  navigator.home();
310  swissCrossVec.push_back(navigator.east());
311  navigator.home();
312  swissCrossVec.push_back(navigator.north());
313  navigator.home();
314  swissCrossVec.push_back(navigator.south());
315  navigator.home();
316 
317  for (unsigned int i = 0; i < swissCrossVec.size(); ++i)
318  {
319 
320  // look for this hit
321  thisHit = recHits_->find(swissCrossVec[i]);
322 
323  // continue if this hit was not found
324  if ((swissCrossVec[i] == DetId(0)) || thisHit == recHits_->end()) continue;
325 
326  // the recHit has to be skipped in the local maximum search if it was found
327  // in the map of channels to be excluded
328  uint32_t rhFlag = thisHit->recoFlag();
329  std::vector<int>::const_iterator vit = std::find(v_chstatus_.begin(), v_chstatus_.end(), rhFlag);
330  if (vit != v_chstatus_.end()) continue;
331 
332  // if this crystal has more energy than the seed then we do
333  // not have a local maxima
334  if (thisHit->energy() > seedEnergy)
335  {
336  maxima = false;
337  break;
338  }
339  }
340 
341  return maxima;
342 
343 }
const EcalRecHitCollection * recHits_
std::vector< EcalRecHit >::const_iterator const_iterator
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
std::vector< int > v_chstatus_
T west() const
move the navigator west
Definition: CaloNavigator.h:59
T south() const
move the navigator south
Definition: CaloNavigator.h:45
T pos() const
get the current position
Definition: CaloNavigator.h:32
const_iterator end() const
T east() const
move the navigator east
Definition: CaloNavigator.h:52
void home() const
move the navigator back to the starting point
Definition: DetId.h:18
iterator find(key_type k)
T north() const
move the navigator north
Definition: CaloNavigator.h:38
void Multi5x5ClusterAlgo::mainSearch ( const EcalRecHitCollection hits,
const CaloSubdetectorGeometry geometry_p,
const CaloSubdetectorTopology topology_p,
const CaloSubdetectorGeometry geometryES_p 
)
private

Definition at line 134 of file Multi5x5ClusterAlgo.cc.

References Multi5x5ClusterAlgo::ProtoBasicCluster::energy(), edm::SortedCollection< T, SORT >::find(), hfClusterShapes_cfi::hits, Multi5x5ClusterAlgo::ProtoBasicCluster::hits(), EcalRecHit::id(), EcalRecHit::kGood, LogTrace, reco::CaloCluster::multi5x5, particleFlowRecHitECAL_cfi::navigator, position, mps_fire::result, Multi5x5ClusterAlgo::ProtoBasicCluster::seed(), and jetUpdater_cfi::sort.

139 {
140 
141  LogTrace("EcalClusters") << "Building clusters............";
142 
143  // Loop over seeds:
144  std::vector<EcalRecHit>::iterator it;
145  for (it = seeds.begin(); it != seeds.end(); it++)
146  {
147 
148  // check if this crystal is able to seed
149  // (event though it is already used)
150  bool usedButCanSeed = false;
151  if (canSeed_s.find(it->id()) != canSeed_s.end()) usedButCanSeed = true;
152 
153  // avoid seeding for anomalous channels
154  if(! it->checkFlag(EcalRecHit::kGood)){ // if rechit is good, no need for further checks
155  if (it->checkFlags( v_chstatus_ )) continue;
156  }
157 
158  // make sure the current seed does not belong to a cluster already.
159  if ((used_s.find(it->id()) != used_s.end()) && (usedButCanSeed == false))
160  {
161  if (it == seeds.begin())
162  {
163  LogTrace("EcalClusters") << "##############################################################" ;
164  LogTrace("EcalClusters") << "DEBUG ALERT: Highest energy seed already belongs to a cluster!";
165  LogTrace("EcalClusters") << "##############################################################";
166 
167  }
168 
169  // seed crystal is used or is used and cannot seed a cluster
170  // so continue to the next seed crystal...
171  continue;
172  }
173 
174  // clear the vector of hits in current cluster
175  current_v.clear();
176 
177  // Create a navigator at the seed and get seed
178  // energy
179  CaloNavigator<DetId> navigator(it->id(), topology_p);
180  DetId seedId = navigator.pos();
181  EcalRecHitCollection::const_iterator seedIt = hits->find(seedId);
182  navigator.setHome(seedId);
183 
184  // Is the seed a local maximum?
185  bool localMaxima = checkMaxima(navigator, hits);
186 
187  if (localMaxima)
188  {
189  // build the 5x5 taking care over which crystals
190  // can seed new clusters and which can't
191  prepareCluster(navigator, hits, geometry_p);
192  }
193 
194  // If some crystals in the current vector then
195  // make them into a cluster
196  if (!current_v.empty())
197  {
198  makeCluster(hits, geometry_p, geometryES_p, seedIt, usedButCanSeed);
199  }
200 
201  } // End loop on seed crystals
202 
203 
206 
207 
208  for(size_t clusNr=0;clusNr<protoClusters_.size();clusNr++){
209  if(!protoClusters_[clusNr].containsSeed()){
210  const EcalRecHit& seedHit =protoClusters_[clusNr].seed();
211  typedef std::vector<std::pair<DetId,int> >::iterator It;
212  std::pair<It,It> result = std::equal_range(whichClusCrysBelongsTo_.begin(),whichClusCrysBelongsTo_.end(),seedHit.id(),PairSortByFirst<DetId,int>());
213 
214  if(result.first!=result.second) protoClusters_[result.first->second].removeHit(seedHit);
215  protoClusters_[clusNr].addSeed();
216 
217  }
218  }
219  }
220 
221 
222 
223  for(size_t clusNr=0;clusNr<protoClusters_.size();clusNr++){
224  const ProtoBasicCluster& protoCluster= protoClusters_[clusNr];
225  Point position;
226  position = posCalculator_.Calculate_Location(protoCluster.hits(), hits, geometry_p, geometryES_p);
227  clusters_v.push_back(reco::BasicCluster(protoCluster.energy(), position, reco::CaloID(detector_), protoCluster.hits(),
228  reco::CaloCluster::multi5x5, protoCluster.seed().id()));
229  }
230 
231  protoClusters_.clear();
232  whichClusCrysBelongsTo_.clear();
233 }
void makeCluster(const EcalRecHitCollection *hits, const CaloSubdetectorGeometry *geometry_p, const CaloSubdetectorGeometry *geometryES_p, const EcalRecHitCollection::const_iterator &seedIt, bool seedOutside)
std::set< DetId > canSeed_s
std::vector< EcalRecHit >::const_iterator const_iterator
std::vector< int > v_chstatus_
void prepareCluster(CaloNavigator< DetId > &navigator, const EcalRecHitCollection *hits, const CaloSubdetectorGeometry *geometry)
std::vector< std::pair< DetId, int > > whichClusCrysBelongsTo_
std::set< DetId > used_s
std::vector< std::pair< DetId, float > > current_v
#define LogTrace(id)
std::vector< ProtoBasicCluster > protoClusters_
std::vector< EcalRecHit > seeds
Definition: DetId.h:18
DetId id() const
get the id
Definition: EcalRecHit.h:77
bool checkMaxima(CaloNavigator< DetId > &navigator, const EcalRecHitCollection *hits)
iterator find(key_type k)
math::XYZPoint Calculate_Location(const HitsAndFractions &iDetIds, const edm::SortedCollection< HitType > *iRecHits, const CaloSubdetectorGeometry *iSubGeom, const CaloSubdetectorGeometry *iESGeom=0)
Definition: PositionCalc.h:68
Structure Point Contains parameters of Gaussian fits to DMRs.
Definition: DMRtrends.cc:55
static int position[264][3]
Definition: ReadPGInfo.cc:509
reco::CaloID::Detectors detector_
The ecal region used.
std::vector< reco::BasicCluster > clusters_v
void Multi5x5ClusterAlgo::makeCluster ( const EcalRecHitCollection hits,
const CaloSubdetectorGeometry geometry_p,
const CaloSubdetectorGeometry geometryES_p,
const EcalRecHitCollection::const_iterator seedIt,
bool  seedOutside 
)
private

Definition at line 235 of file Multi5x5ClusterAlgo.cc.

References reco::CaloID::DET_ECAL_BARREL, reco::CaloID::DET_ECAL_ENDCAP, EcalBarrel, EcalRecHit::energy(), edm::SortedCollection< T, SORT >::find(), LogTrace, and position.

240 {
241 
242  double energy = 0;
243  //double chi2 = 0;
244  reco::CaloID caloID;
245  Point position;
246  position = posCalculator_.Calculate_Location(current_v, hits, geometry, geometryES);
247 
248  std::vector<std::pair<DetId, float> >::iterator it;
249  for (it = current_v.begin(); it != current_v.end(); it++)
250  {
251  EcalRecHitCollection::const_iterator itt = hits->find( (*it).first );
252  EcalRecHit hit_p = *itt;
253  energy += hit_p.energy();
254  //chi2 += 0;
255  if ( (*it).first.subdetId() == EcalBarrel ) {
257  } else {
259  }
260 
261  }
262  //chi2 /= energy;
263 
264  LogTrace("EcalClusters") << "******** NEW CLUSTER ********";
265  LogTrace("EcalClusters") << "No. of crystals = " << current_v.size();
266  LogTrace("EcalClusters") << " Energy = " << energy ;
267  LogTrace("EcalClusters") << " Phi = " << position.phi();
268  LogTrace("EcalClusters") << " Eta " << position.eta();
269  LogTrace("EcalClusters") << "*****************************";
270 
271 
272  // to be a valid cluster the cluster energy
273  // must be at least the seed energy
274  double seedEnergy = seedIt->energy();
275  if ((seedOutside && energy>=0) || (!seedOutside && energy >= seedEnergy))
276  {
277  if(reassignSeedCrysToClusterItSeeds_){ //if we're not doing this, we dont need this info so lets not bother filling it
278  for(size_t hitNr=0;hitNr<current_v.size();hitNr++) whichClusCrysBelongsTo_.push_back(std::pair<DetId,int>(current_v[hitNr].first,protoClusters_.size()));
279  }
280  protoClusters_.push_back(ProtoBasicCluster(energy,*seedIt,current_v));
281 
282  // clusters_v.push_back(reco::BasicCluster(energy, position, reco::CaloID(detector_), current_v, reco::CaloCluster::multi5x5, seedIt->id()));
283 
284  // if no valid cluster was built,
285  // then free up these crystals to be used in the next...
286  } else {
287  std::vector<std::pair<DetId, float> >::iterator iter;
288  for (iter = current_v.begin(); iter != current_v.end(); iter++)
289  {
290  used_s.erase(iter->first);
291  } //for(iter)
292  } //else
293 
294 }
std::vector< EcalRecHit >::const_iterator const_iterator
std::vector< std::pair< DetId, int > > whichClusCrysBelongsTo_
std::set< DetId > used_s
float energy() const
Definition: EcalRecHit.h:68
std::vector< std::pair< DetId, float > > current_v
#define LogTrace(id)
std::vector< ProtoBasicCluster > protoClusters_
iterator find(key_type k)
math::XYZPoint Calculate_Location(const HitsAndFractions &iDetIds, const edm::SortedCollection< HitType > *iRecHits, const CaloSubdetectorGeometry *iSubGeom, const CaloSubdetectorGeometry *iESGeom=0)
Definition: PositionCalc.h:68
Structure Point Contains parameters of Gaussian fits to DMRs.
Definition: DMRtrends.cc:55
static int position[264][3]
Definition: ReadPGInfo.cc:509
std::vector< reco::BasicCluster > Multi5x5ClusterAlgo::makeClusters ( const EcalRecHitCollection hits,
const CaloSubdetectorGeometry geometry,
const CaloSubdetectorTopology topology_p,
const CaloSubdetectorGeometry geometryES_p,
reco::CaloID::Detectors  detector,
bool  regional = false,
const std::vector< RectangularEtaPhiRegion > &  regions = std::vector<RectangularEtaPhiRegion>() 
)

Definition at line 46 of file Multi5x5ClusterAlgo.cc.

References edm::SortedCollection< T, SORT >::begin(), reco::CaloID::DET_ECAL_BARREL, reco::CaloID::DET_ECAL_ENDCAP, reco::CaloID::DET_NONE, edm::SortedCollection< T, SORT >::end(), ET, CaloSubdetectorGeometry::getGeometry(), hfClusterShapes_cfi::hits, isClusterEtLess(), LogTrace, jetUpdater_cfi::sort, AlCaHLTBitMon_QueryRunRegistry::string, electronIdCutBased_cfi::threshold, x, and y.

Referenced by Multi5x5ClusterProducer::clusterizeECALPart(), EgammaHLTMulti5x5ClusterProducer::clusterizeECALPart(), and ~Multi5x5ClusterAlgo().

55 {
56  seeds.clear();
57  used_s.clear();
58  canSeed_s.clear();
59  clusters_v.clear();
60 
61  recHits_ = hits;
62 
63  double threshold = 0;
64  std::string ecalPart_string;
67  {
69  threshold = ecalEndcapSeedThreshold;
70  ecalPart_string = "EndCap";
71  }
73  {
75  threshold = ecalBarrelSeedThreshold;
76  ecalPart_string = "Barrel";
77  }
78 
79  LogTrace("EcalClusters") << "-------------------------------------------------------------";
80  LogTrace("EcalClusters") << "Island algorithm invoked for ECAL" << ecalPart_string ;
81  LogTrace("EcalClusters") << "Looking for seeds, energy threshold used = " << threshold << " GeV";
82 
83 
84  int nregions=0;
85  if(regional) nregions=regions.size();
86 
87  if(!regional || nregions) {
88 
90  for(it = hits->begin(); it != hits->end(); it++)
91  {
92  double energy = it->energy();
93  if (energy < threshold) continue; // need to check to see if this line is useful!
94 
95  auto thisCell = geometry_p->getGeometry(it->id());
96  // Require that RecHit is within clustering region in case
97  // of regional reconstruction
98  bool withinRegion = false;
99  if (regional) {
100  std::vector<RectangularEtaPhiRegion>::const_iterator region;
101  for (region=regions.begin(); region!=regions.end(); region++) {
102  if (region->inRegion(thisCell->etaPos(),thisCell->phiPos())) {
103  withinRegion = true;
104  break;
105  }
106  }
107  }
108 
109  if (!regional || withinRegion) {
110  float ET = it->energy() * thisCell->getPosition().basicVector().unit().perp();
111  if (ET > threshold) seeds.push_back(*it);
112  }
113  }
114 
115  }
116 
117  sort(seeds.begin(), seeds.end(), [](auto const& x, auto const& y){return (x.energy() > y.energy());});
118 
119  LogTrace("EcalClusters") << "Total number of seeds found in event = " << seeds.size();
120 
121 
122  mainSearch(hits, geometry_p, topology_p, geometryES_p);
123  sort(clusters_v.rbegin(), clusters_v.rend(), isClusterEtLess);
124 
125  LogTrace("EcalClusters") << "---------- end of main search. clusters have been sorted ----";
126 
127 
128  return clusters_v;
129 
130 }
std::set< DetId > canSeed_s
const EcalRecHitCollection * recHits_
std::vector< EcalRecHit >::const_iterator const_iterator
void mainSearch(const EcalRecHitCollection *hits, const CaloSubdetectorGeometry *geometry_p, const CaloSubdetectorTopology *topology_p, const CaloSubdetectorGeometry *geometryES_p)
std::set< DetId > used_s
bool isClusterEtLess(const reco::CaloCluster &x, const reco::CaloCluster &y)
Definition: ClusterEtLess.h:7
#define LogTrace(id)
const_iterator end() const
std::vector< EcalRecHit > seeds
reco::CaloID::Detectors detector_
The ecal region used.
std::vector< reco::BasicCluster > clusters_v
#define ET
const_iterator begin() const
void Multi5x5ClusterAlgo::prepareCluster ( CaloNavigator< DetId > &  navigator,
const EcalRecHitCollection hits,
const CaloSubdetectorGeometry geometry 
)
private

Definition at line 345 of file Multi5x5ClusterAlgo.cc.

References funct::abs(), PVValHelper::dx, PVValHelper::dy, CaloNavigator< T, TOPO >::home(), and CaloNavigator< T, TOPO >::offsetBy().

348 {
349 
350  DetId thisDet;
351  std::set<DetId>::iterator setItr;
352 
353  // now add the 5x5 taking care to mark the edges
354  // as able to seed and where overlapping in the central
355  // region with crystals that were previously able to seed
356  // change their status so they are not able to seed
357  //std::cout << std::endl;
358  for (int dx = -2; dx < 3; ++dx)
359  {
360  for (int dy = -2; dy < 3; ++ dy)
361  {
362 
363  // navigate in free steps forming
364  // a full 5x5
365  thisDet = navigator.offsetBy(dx, dy);
366  navigator.home();
367 
368  // add the current crystal
369  //std::cout << "adding " << dx << ", " << dy << std::endl;
370  addCrystal(thisDet);
371 
372  // now consider if we are in an edge (outer 16)
373  // or central (inner 9) region
374  if ((abs(dx) > 1) || (abs(dy) > 1))
375  {
376  // this is an "edge" so should be allowed to seed
377  // provided it is not already used
378  //std::cout << " setting can seed" << std::endl;
379  canSeed_s.insert(thisDet);
380  } // end if "edge"
381  else
382  {
383  // or else we are in the central 3x3
384  // and must remove any of these crystals from the canSeed set
385  setItr = canSeed_s.find(thisDet);
386  if (setItr != canSeed_s.end())
387  {
388  //std::cout << " unsetting can seed" << std::endl;
389  canSeed_s.erase(setItr);
390  }
391  } // end if "centre"
392 
393 
394  } // end loop on dy
395 
396  } // end loop on dx
397 
398  //std::cout << "*** " << std::endl;
399  //std::cout << " current_v contains " << current_v.size() << std::endl;
400  //std::cout << "*** " << std::endl;
401 }
std::set< DetId > canSeed_s
T offsetBy(int deltaX, int deltaY) const
Free movement of arbitray steps.
Definition: CaloNavigator.h:80
void addCrystal(const DetId &det)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void home() const
move the navigator back to the starting point
Definition: DetId.h:18

Member Data Documentation

std::set<DetId> Multi5x5ClusterAlgo::canSeed_s
private

Definition at line 104 of file Multi5x5ClusterAlgo.h.

std::vector<reco::BasicCluster> Multi5x5ClusterAlgo::clusters_v
private

Definition at line 112 of file Multi5x5ClusterAlgo.h.

std::vector<std::pair<DetId, float> > Multi5x5ClusterAlgo::current_v
private

Definition at line 109 of file Multi5x5ClusterAlgo.h.

reco::CaloID::Detectors Multi5x5ClusterAlgo::detector_
private

The ecal region used.

Definition at line 88 of file Multi5x5ClusterAlgo.h.

double Multi5x5ClusterAlgo::ecalBarrelSeedThreshold
private

Definition at line 91 of file Multi5x5ClusterAlgo.h.

double Multi5x5ClusterAlgo::ecalEndcapSeedThreshold
private

Definition at line 92 of file Multi5x5ClusterAlgo.h.

PositionCalc Multi5x5ClusterAlgo::posCalculator_
private

Definition at line 85 of file Multi5x5ClusterAlgo.h.

Referenced by Multi5x5ClusterAlgo().

std::vector<ProtoBasicCluster> Multi5x5ClusterAlgo::protoClusters_
private

Definition at line 113 of file Multi5x5ClusterAlgo.h.

bool Multi5x5ClusterAlgo::reassignSeedCrysToClusterItSeeds_
private

Definition at line 117 of file Multi5x5ClusterAlgo.h.

const EcalRecHitCollection* Multi5x5ClusterAlgo::recHits_
private

Definition at line 95 of file Multi5x5ClusterAlgo.h.

std::vector<EcalRecHit> Multi5x5ClusterAlgo::seeds
private

Definition at line 98 of file Multi5x5ClusterAlgo.h.

std::set<DetId> Multi5x5ClusterAlgo::used_s
private

Definition at line 103 of file Multi5x5ClusterAlgo.h.

std::vector<int> Multi5x5ClusterAlgo::v_chstatus_
private

Definition at line 115 of file Multi5x5ClusterAlgo.h.

Referenced by Multi5x5ClusterAlgo().

std::vector<std::pair<DetId,int> > Multi5x5ClusterAlgo::whichClusCrysBelongsTo_
private

Definition at line 100 of file Multi5x5ClusterAlgo.h.