CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Classes | Public Types | Public Member Functions | Private Member Functions | Private Attributes
PFECALSuperClusterAlgo Class Reference

#include <PFECALSuperClusterAlgo.h>

Classes

class  CalibratedPFCluster
 

Public Types

typedef std::shared_ptr
< CalibratedPFCluster
CalibratedClusterPtr
 
typedef std::vector
< CalibratedClusterPtr
CalibratedClusterPtrVector
 
enum  clustering_type { kBOX =1, kMustache =2 }
 

Public Member Functions

std::auto_ptr
< reco::SuperClusterCollection
getEBOutputSCCollection ()
 
std::auto_ptr
< reco::SuperClusterCollection
getEEOutputSCCollection ()
 
void loadAndSortPFClusters (const edm::View< reco::PFCluster > &ecalclusters, const edm::View< reco::PFCluster > &psclusters)
 
 PFECALSuperClusterAlgo ()
 constructor More...
 
void run ()
 
void setClusteringType (clustering_type thetype)
 
void setCrackCorrections (bool applyCrackCorrections)
 
void setEtawidthSuperClusterBarrel (double etawidth)
 
void setEtawidthSuperClusterEndcap (double etawidth)
 
void setMajorityFraction (const double f)
 
void setPFClusterCalibration (const std::shared_ptr< PFEnergyCalibration > &)
 
void setPhiwidthSuperClusterBarrel (double phiwidth)
 
void setPhiwidthSuperClusterEndcap (double phiwidth)
 
void setSatelliteMerging (const bool doit)
 
void setSatelliteThreshold (const double t)
 
void setThreshPFClusterBarrel (double thresh)
 
void setThreshPFClusterEndcap (double thresh)
 
void setThreshPFClusterES (double thresh)
 
void setThreshPFClusterSeedBarrel (double thresh)
 
void setThreshPFClusterSeedEndcap (double thresh)
 
void setUseDynamicDPhi (bool useit)
 
void setUsePS (bool useit)
 
void setVerbosityLevel (bool verbose)
 

Private Member Functions

void buildAllSuperClusters (CalibratedClusterPtrVector &, double seedthresh)
 
void buildSuperCluster (CalibratedClusterPtr &, CalibratedClusterPtrVector &)
 

Private Attributes

CalibratedClusterPtrVector _clustersEB
 
CalibratedClusterPtrVector _clustersEE
 
clustering_type _clustype
 
std::shared_ptr
< PFEnergyCalibration
_pfEnergyCalibration
 
std::unordered_map< edm::Ptr
< reco::PFCluster >
, edm::PtrVector
< reco::PFCluster > > 
_psclustersforee
 
bool _useDynamicDPhi
 
bool applyCrackCorrections_
 
bool doSatelliteClusterMerge_
 
double etawidthSuperCluster_
 
double etawidthSuperClusterBarrel_
 
double etawidthSuperClusterEndcap_
 
double fractionForMajority_
 
double phiwidthSuperCluster_
 
double phiwidthSuperClusterBarrel_
 
double phiwidthSuperClusterEndcap_
 
double satelliteThreshold_
 
std::auto_ptr
< reco::SuperClusterCollection
superClustersEB_
 
std::auto_ptr
< reco::SuperClusterCollection
superClustersEE_
 
double threshPFCluster_
 
double threshPFClusterBarrel_
 
double threshPFClusterEndcap_
 
double threshPFClusterES_
 
double threshPFClusterSeed_
 
double threshPFClusterSeedBarrel_
 
double threshPFClusterSeedEndcap_
 
bool usePS
 
bool verbose_
 

Detailed Description

Definition at line 47 of file PFECALSuperClusterAlgo.h.

Member Typedef Documentation

Definition at line 70 of file PFECALSuperClusterAlgo.h.

Definition at line 71 of file PFECALSuperClusterAlgo.h.

Member Enumeration Documentation

Constructor & Destructor Documentation

PFECALSuperClusterAlgo::PFECALSuperClusterAlgo ( )

constructor

Definition at line 172 of file PFECALSuperClusterAlgo.cc.

172 { }

Member Function Documentation

void PFECALSuperClusterAlgo::buildAllSuperClusters ( CalibratedClusterPtrVector ,
double  seedthresh 
)
private

Definition at line 259 of file PFECALSuperClusterAlgo.cc.

260  {
261  IsASeed seedable(seedthresh);
262  // make sure only seeds appear at the front of the list of clusters
263  std::stable_partition(clusters.begin(),clusters.end(),seedable);
264  // in each iteration we are working on a list that is already sorted
265  // in the cluster energy and remains so through each iteration
266  // NB: since clusters is sorted in loadClusters any_of has O(1)
267  // timing until you run out of seeds!
268  while( std::any_of(clusters.cbegin(), clusters.cend(), seedable) ) {
269  buildSuperCluster(clusters.front(),clusters);
270  }
271 }
void buildSuperCluster(CalibratedClusterPtr &, CalibratedClusterPtrVector &)
void PFECALSuperClusterAlgo::buildSuperCluster ( CalibratedClusterPtr ,
CalibratedClusterPtrVector  
)
private

Definition at line 274 of file PFECALSuperClusterAlgo.cc.

References reco::SuperCluster::addCluster(), reco::CaloCluster::addHitAndFraction(), reco::SuperCluster::addPreshowerCluster(), PFLayer::ECAL_BARREL, PFLayer::ECAL_ENDCAP, edm::hlt::Exception, spr::find(), edm::Ptr< T >::get(), PFClusterWidthAlgo::pflowEtaWidth(), PFClusterWidthAlgo::pflowPhiWidth(), reco::SuperCluster::preshowerClustersBegin(), reco::SuperCluster::preshowerClustersEnd(), PFLayer::PS1, PFLayer::PS2, reco::SuperCluster::rawEnergy(), reco::SuperCluster::setEtaWidth(), reco::SuperCluster::setPhiWidth(), reco::SuperCluster::setPreshowerEnergy(), reco::SuperCluster::setPreshowerEnergyPlane1(), reco::SuperCluster::setPreshowerEnergyPlane2(), and reco::SuperCluster::setSeed().

275  {
276  IsClustered IsClusteredWithSeed(seed,_clustype,_useDynamicDPhi);
277  IsLinkedByRecHit MatchesSeedByRecHit(seed,satelliteThreshold_,
278  fractionForMajority_,0.1,0.2);
279  bool isEE = false;
280  SumPSEnergy sumps1(PFLayer::PS1), sumps2(PFLayer::PS2);
281  switch( seed->the_ptr()->layer() ) {
283  IsClusteredWithSeed.phiwidthSuperCluster_ = phiwidthSuperClusterBarrel_;
284  IsClusteredWithSeed.etawidthSuperCluster_ = etawidthSuperClusterBarrel_;
285  edm::LogInfo("PFClustering") << "Building SC number "
286  << superClustersEB_->size() + 1
287  << " in the ECAL barrel!";
288  break;
289  case PFLayer::ECAL_ENDCAP:
290  IsClusteredWithSeed.phiwidthSuperCluster_ = phiwidthSuperClusterEndcap_;
291  IsClusteredWithSeed.etawidthSuperCluster_ = etawidthSuperClusterEndcap_;
292  edm::LogInfo("PFClustering") << "Building SC number "
293  << superClustersEE_->size() + 1
294  << " in the ECAL endcap!" << std::endl;
295  isEE = true;
296  break;
297  default:
298  break;
299  }
300 
301  // this function shuffles the list of clusters into a list
302  // where all clustered sub-clusters are at the front
303  // and returns a pointer to the first unclustered cluster.
304  // The relative ordering of clusters is preserved
305  // (i.e. both resulting sub-lists are sorted by energy).
306  auto not_clustered = std::stable_partition(clusters.begin(),clusters.end(),
307  IsClusteredWithSeed);
308  // satellite cluster merging
309  // it was found that large clusters can split!
310  if( doSatelliteClusterMerge_ ) {
311  not_clustered = std::stable_partition(not_clustered,clusters.end(),
312  MatchesSeedByRecHit);
313  }
314 
315  if(verbose_) {
316  edm::LogInfo("PFClustering") << "Dumping cluster detail";
317  edm::LogVerbatim("PFClustering")
318  << "\tPassed seed: e = " << seed->energy_nocalib()
319  << " eta = " << seed->eta() << " phi = " << seed->phi()
320  << std::endl;
321  for( auto clus = clusters.cbegin(); clus != not_clustered; ++clus ) {
322  edm::LogVerbatim("PFClustering")
323  << "\t\tClustered cluster: e = " << (*clus)->energy_nocalib()
324  << " eta = " << (*clus)->eta() << " phi = " << (*clus)->phi()
325  << std::endl;
326  }
327  for( auto clus = not_clustered; clus != clusters.end(); ++clus ) {
328  edm::LogVerbatim("PFClustering")
329  << "\tNon-Clustered cluster: e = " << (*clus)->energy_nocalib()
330  << " eta = " << (*clus)->eta() << " phi = " << (*clus)->phi()
331  << std::endl;
332  }
333  }
334  // move the clustered clusters out of available cluster list
335  // and into a temporary vector for building the SC
336  CalibratedClusterPtrVector clustered(clusters.begin(),not_clustered);
337  clusters.erase(clusters.begin(),not_clustered);
338  // need the vector of raw pointers for a PF width class
339  std::vector<const reco::PFCluster*> bare_ptrs;
340  // calculate necessary parameters and build the SC
341  double posX(0), posY(0), posZ(0),
342  rawSCEnergy(0), corrSCEnergy(0), clusterCorrEE(0),
343  PS1_clus_sum(0), PS2_clus_sum(0);
344  for( auto& clus : clustered ) {
345  bare_ptrs.push_back(clus->the_ptr().get());
346 
347  const double cluseraw = clus->energy_nocalib();
348  const math::XYZPoint& cluspos = clus->the_ptr()->position();
349  posX += cluseraw * cluspos.X();
350  posY += cluseraw * cluspos.Y();
351  posZ += cluseraw * cluspos.Z();
352  // update EE calibrated super cluster energies
353  if( isEE ) {
354  const auto& psclusters = _psclustersforee[clus->the_ptr()];
355  PS1_clus_sum = std::accumulate(psclusters.begin(),psclusters.end(),
356  0.0,sumps1);
357  PS2_clus_sum = std::accumulate(psclusters.begin(),psclusters.end(),
358  0.0,sumps2);
359  clusterCorrEE =
360  _pfEnergyCalibration->energyEm(*(clus->the_ptr()),
361  PS1_clus_sum,PS2_clus_sum,
363  clus->resetCalibratedEnergy(clusterCorrEE);
364  }
365 
366  rawSCEnergy += cluseraw;
367  corrSCEnergy += clus->energy();
368  }
369  posX /= rawSCEnergy;
370  posY /= rawSCEnergy;
371  posZ /= rawSCEnergy;
372 
373  // now build the supercluster
374  reco::SuperCluster new_sc(corrSCEnergy,math::XYZPoint(posX,posY,posZ));
375  double ps1_energy(0.0), ps2_energy(0.0), ps_energy(0.0);
376  new_sc.setSeed(clustered.front()->the_ptr());
377  for( const auto& clus : clustered ) {
378  new_sc.addCluster(clus->the_ptr());
379  auto& hits_and_fractions = clus->the_ptr()->hitsAndFractions();
380  for( auto& hit_and_fraction : hits_and_fractions ) {
381  new_sc.addHitAndFraction(hit_and_fraction.first,hit_and_fraction.second);
382  }
383  const auto& cluspsassociation = _psclustersforee[clus->the_ptr()];
384  // EE rechits should be uniquely matched to sets of pre-shower
385  // clusters at this point, so we throw an exception if otherwise
386  for( const auto& psclus : cluspsassociation ) {
387  auto found_pscluster = std::find(new_sc.preshowerClustersBegin(),
388  new_sc.preshowerClustersEnd(),
389  reco::CaloClusterPtr(psclus));
390  if( found_pscluster == new_sc.preshowerClustersEnd() ) {
391  const double psenergy = psclus->energy();
392  new_sc.addPreshowerCluster(psclus);
393  ps1_energy += (PFLayer::PS1 == psclus->layer())*psenergy;
394  ps2_energy += (PFLayer::PS2 == psclus->layer())*psenergy;
395  ps_energy += psenergy;
396  } else {
397  throw cms::Exception("PFECALSuperClusterAlgo::buildSuperCluster")
398  << "Found a PS cluster matched to more than one EE cluster!"
399  << std::endl << std::hex << psclus.get() << " == "
400  << found_pscluster->get() << std::dec << std::endl;
401  }
402  }
403  }
404  new_sc.setPreshowerEnergy(ps_energy);
405  new_sc.setPreshowerEnergyPlane1(ps1_energy);
406  new_sc.setPreshowerEnergyPlane2(ps2_energy);
407 
408  // calculate linearly weighted cluster widths
409  PFClusterWidthAlgo pfwidth(bare_ptrs);
410  new_sc.setEtaWidth(pfwidth.pflowEtaWidth());
411  new_sc.setPhiWidth(pfwidth.pflowPhiWidth());
412 
413  // cache the value of the raw energy
414  new_sc.rawEnergy();
415 
416  // save the super cluster to the appropriate list
417  switch( seed->the_ptr()->layer() ) {
419  superClustersEB_->push_back(new_sc);
420  break;
421  case PFLayer::ECAL_ENDCAP:
422  superClustersEE_->push_back(new_sc);
423  break;
424  default:
425  break;
426  }
427 }
std::auto_ptr< reco::SuperClusterCollection > superClustersEB_
std::unordered_map< edm::Ptr< reco::PFCluster >, edm::PtrVector< reco::PFCluster > > _psclustersforee
edm::Ptr< CaloCluster > CaloClusterPtr
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
std::vector< CalibratedClusterPtr > CalibratedClusterPtrVector
std::shared_ptr< PFEnergyCalibration > _pfEnergyCalibration
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:13
std::auto_ptr< reco::SuperClusterCollection > superClustersEE_
std::auto_ptr<reco::SuperClusterCollection> PFECALSuperClusterAlgo::getEBOutputSCCollection ( )
inline

Definition at line 107 of file PFECALSuperClusterAlgo.h.

References superClustersEB_.

107 { return superClustersEB_; }
std::auto_ptr< reco::SuperClusterCollection > superClustersEB_
std::auto_ptr<reco::SuperClusterCollection> PFECALSuperClusterAlgo::getEEOutputSCCollection ( )
inline

Definition at line 109 of file PFECALSuperClusterAlgo.h.

References superClustersEE_.

109 { return superClustersEE_; }
std::auto_ptr< reco::SuperClusterCollection > superClustersEE_
void PFECALSuperClusterAlgo::loadAndSortPFClusters ( const edm::View< reco::PFCluster > &  ecalclusters,
const edm::View< reco::PFCluster > &  psclusters 
)

Definition at line 180 of file PFECALSuperClusterAlgo.cc.

References PFLayer::ECAL_BARREL, PFLayer::ECAL_ENDCAP, edm::Ptr< T >::isNonnull(), LogDebug, PFLayer::PS1, PFLayer::PS2, and python.multivaluedict::sort().

181  {
182  // reset the system for running
184  _clustersEB.clear();
186  _clustersEE.clear();
187  _psclustersforee.clear();
188 
189  auto clusterPtrs = clusters.ptrVector();
190  //Select PF clusters available for the clustering
191  for ( auto& cluster : clusterPtrs ){
192  LogDebug("PFClustering")
193  << "Loading PFCluster i="<<cluster.key()
194  <<" energy="<<cluster->energy()<<std::endl;
195 
196  double Ecorr = _pfEnergyCalibration->energyEm(*cluster,
197  0.0,0.0,
199  CalibratedClusterPtr calib_cluster(new CalibratedPFCluster(cluster,Ecorr));
200  switch( cluster->layer() ) {
202  if( calib_cluster->energy() > threshPFClusterBarrel_ ) {
203  _clustersEB.push_back(calib_cluster);
204  }
205  break;
207  if( calib_cluster->energy() > threshPFClusterEndcap_ ) {
208  _clustersEE.push_back(calib_cluster);
209  _psclustersforee.emplace(calib_cluster->the_ptr(),
211  }
212  break;
213  default:
214  break;
215  }
216  }
217  // make the association map of ECAL clusters to preshower clusters
218  edm::PtrVector<reco::PFCluster> clusterPtrsPS = psclusters.ptrVector();
219  double dist = -1.0, min_dist = -1.0;
220  // match PS clusters to EE clusters, minimum distance to EE is ensured
221  // since the inner loop is over the EE clusters
222  for( const auto& psclus : clusterPtrsPS ) {
223  if( psclus->energy() < threshPFClusterES_ ) continue;
224  switch( psclus->layer() ) { // just in case this isn't the ES...
225  case PFLayer::PS1:
226  case PFLayer::PS2:
227  break;
228  default:
229  continue;
230  }
232  dist = min_dist = -1.0; // reset
233  for( const auto& eeclus : _clustersEE ) {
234  dist = testPreshowerDistance(eeclus->the_ptr(),psclus);
235  if( dist == -1.0 || (min_dist != -1.0 && dist > min_dist) ) continue;
236  if( dist < min_dist || min_dist == -1.0 ) {
237  eematch = eeclus->the_ptr();
238  min_dist = dist;
239  }
240  } // loop on EE clusters
241  if( eematch.isNonnull() ) _psclustersforee[eematch].push_back(psclus);
242  } // loop on PS clusters
243 
244  // sort full cluster collections by their calibrated energy
245  // this will put all the seeds first by construction
246  GreaterByEt greater;
247  std::sort(_clustersEB.begin(), _clustersEB.end(), greater);
248  std::sort(_clustersEE.begin(), _clustersEE.end(), greater);
249 }
#define LogDebug(id)
std::auto_ptr< reco::SuperClusterCollection > superClustersEB_
std::unordered_map< edm::Ptr< reco::PFCluster >, edm::PtrVector< reco::PFCluster > > _psclustersforee
std::shared_ptr< CalibratedPFCluster > CalibratedClusterPtr
bool isNonnull() const
Checks for non-null.
Definition: Ptr.h:152
std::shared_ptr< PFEnergyCalibration > _pfEnergyCalibration
std::vector< SuperCluster > SuperClusterCollection
collection of SuperCluser objectr
CalibratedClusterPtrVector _clustersEE
std::auto_ptr< reco::SuperClusterCollection > superClustersEE_
PtrVector< T > const & ptrVector() const
Definition: View.h:125
CalibratedClusterPtrVector _clustersEB
void PFECALSuperClusterAlgo::run ( void  )

Definition at line 251 of file PFECALSuperClusterAlgo.cc.

251  {
252  // clusterize the EB
254  // clusterize the EE
256 }
void buildAllSuperClusters(CalibratedClusterPtrVector &, double seedthresh)
CalibratedClusterPtrVector _clustersEE
CalibratedClusterPtrVector _clustersEB
void PFECALSuperClusterAlgo::setClusteringType ( clustering_type  thetype)
inline

Definition at line 79 of file PFECALSuperClusterAlgo.h.

References _clustype.

79 { _clustype = thetype; }
void PFECALSuperClusterAlgo::setCrackCorrections ( bool  applyCrackCorrections)
inline

Definition at line 104 of file PFECALSuperClusterAlgo.h.

References applyCrackCorrections_.

104 { applyCrackCorrections_ = applyCrackCorrections;}
void PFECALSuperClusterAlgo::setEtawidthSuperClusterBarrel ( double  etawidth)
inline

Definition at line 89 of file PFECALSuperClusterAlgo.h.

References etawidthSuperClusterBarrel_.

void PFECALSuperClusterAlgo::setEtawidthSuperClusterEndcap ( double  etawidth)
inline

Definition at line 91 of file PFECALSuperClusterAlgo.h.

References etawidthSuperClusterEndcap_.

void PFECALSuperClusterAlgo::setMajorityFraction ( const double  f)
inline

Definition at line 100 of file PFECALSuperClusterAlgo.h.

References f, and fractionForMajority_.

100 { fractionForMajority_ = f; }
double f[11][100]
void PFECALSuperClusterAlgo::setPFClusterCalibration ( const std::shared_ptr< PFEnergyCalibration > &  calib)

Definition at line 175 of file PFECALSuperClusterAlgo.cc.

References calib.

175  {
177 }
MVATrainerComputer * calib
Definition: MVATrainer.cc:64
std::shared_ptr< PFEnergyCalibration > _pfEnergyCalibration
void PFECALSuperClusterAlgo::setPhiwidthSuperClusterBarrel ( double  phiwidth)
inline

Definition at line 88 of file PFECALSuperClusterAlgo.h.

References phiwidthSuperClusterBarrel_.

void PFECALSuperClusterAlgo::setPhiwidthSuperClusterEndcap ( double  phiwidth)
inline

Definition at line 90 of file PFECALSuperClusterAlgo.h.

References phiwidthSuperClusterEndcap_.

void PFECALSuperClusterAlgo::setSatelliteMerging ( const bool  doit)
inline

Definition at line 98 of file PFECALSuperClusterAlgo.h.

References doSatelliteClusterMerge_.

void PFECALSuperClusterAlgo::setSatelliteThreshold ( const double  t)
inline

Definition at line 99 of file PFECALSuperClusterAlgo.h.

References satelliteThreshold_, and lumiQTWidget::t.

void PFECALSuperClusterAlgo::setThreshPFClusterBarrel ( double  thresh)
inline
void PFECALSuperClusterAlgo::setThreshPFClusterEndcap ( double  thresh)
inline
void PFECALSuperClusterAlgo::setThreshPFClusterES ( double  thresh)
inline
void PFECALSuperClusterAlgo::setThreshPFClusterSeedBarrel ( double  thresh)
inline
void PFECALSuperClusterAlgo::setThreshPFClusterSeedEndcap ( double  thresh)
inline
void PFECALSuperClusterAlgo::setUseDynamicDPhi ( bool  useit)
inline

Definition at line 81 of file PFECALSuperClusterAlgo.h.

References _useDynamicDPhi.

void PFECALSuperClusterAlgo::setUsePS ( bool  useit)
inline

Definition at line 92 of file PFECALSuperClusterAlgo.h.

References usePS.

92 { usePS = useit; }
void PFECALSuperClusterAlgo::setVerbosityLevel ( bool  verbose)
inline

Member Data Documentation

CalibratedClusterPtrVector PFECALSuperClusterAlgo::_clustersEB
private

Definition at line 118 of file PFECALSuperClusterAlgo.h.

CalibratedClusterPtrVector PFECALSuperClusterAlgo::_clustersEE
private

Definition at line 119 of file PFECALSuperClusterAlgo.h.

clustering_type PFECALSuperClusterAlgo::_clustype
private

Definition at line 125 of file PFECALSuperClusterAlgo.h.

Referenced by setClusteringType().

std::shared_ptr<PFEnergyCalibration> PFECALSuperClusterAlgo::_pfEnergyCalibration
private

Definition at line 124 of file PFECALSuperClusterAlgo.h.

std::unordered_map<edm::Ptr<reco::PFCluster>, edm::PtrVector<reco::PFCluster> > PFECALSuperClusterAlgo::_psclustersforee
private

Definition at line 121 of file PFECALSuperClusterAlgo.h.

bool PFECALSuperClusterAlgo::_useDynamicDPhi
private

Definition at line 152 of file PFECALSuperClusterAlgo.h.

Referenced by setUseDynamicDPhi().

bool PFECALSuperClusterAlgo::applyCrackCorrections_
private

Definition at line 154 of file PFECALSuperClusterAlgo.h.

Referenced by setCrackCorrections().

bool PFECALSuperClusterAlgo::doSatelliteClusterMerge_
private

Definition at line 149 of file PFECALSuperClusterAlgo.h.

Referenced by setSatelliteMerging().

double PFECALSuperClusterAlgo::etawidthSuperCluster_
private

Definition at line 135 of file PFECALSuperClusterAlgo.h.

double PFECALSuperClusterAlgo::etawidthSuperClusterBarrel_
private

Definition at line 145 of file PFECALSuperClusterAlgo.h.

Referenced by setEtawidthSuperClusterBarrel().

double PFECALSuperClusterAlgo::etawidthSuperClusterEndcap_
private

Definition at line 147 of file PFECALSuperClusterAlgo.h.

Referenced by setEtawidthSuperClusterEndcap().

double PFECALSuperClusterAlgo::fractionForMajority_
private

Definition at line 150 of file PFECALSuperClusterAlgo.h.

Referenced by setMajorityFraction().

double PFECALSuperClusterAlgo::phiwidthSuperCluster_
private

Definition at line 136 of file PFECALSuperClusterAlgo.h.

double PFECALSuperClusterAlgo::phiwidthSuperClusterBarrel_
private

Definition at line 144 of file PFECALSuperClusterAlgo.h.

Referenced by setPhiwidthSuperClusterBarrel().

double PFECALSuperClusterAlgo::phiwidthSuperClusterEndcap_
private

Definition at line 146 of file PFECALSuperClusterAlgo.h.

Referenced by setPhiwidthSuperClusterEndcap().

double PFECALSuperClusterAlgo::satelliteThreshold_
private

Definition at line 150 of file PFECALSuperClusterAlgo.h.

Referenced by setSatelliteThreshold().

std::auto_ptr<reco::SuperClusterCollection> PFECALSuperClusterAlgo::superClustersEB_
private

Definition at line 122 of file PFECALSuperClusterAlgo.h.

Referenced by getEBOutputSCCollection().

std::auto_ptr<reco::SuperClusterCollection> PFECALSuperClusterAlgo::superClustersEE_
private

Definition at line 123 of file PFECALSuperClusterAlgo.h.

Referenced by getEEOutputSCCollection().

double PFECALSuperClusterAlgo::threshPFCluster_
private

Definition at line 134 of file PFECALSuperClusterAlgo.h.

double PFECALSuperClusterAlgo::threshPFClusterBarrel_
private

Definition at line 139 of file PFECALSuperClusterAlgo.h.

Referenced by setThreshPFClusterBarrel().

double PFECALSuperClusterAlgo::threshPFClusterEndcap_
private

Definition at line 141 of file PFECALSuperClusterAlgo.h.

Referenced by setThreshPFClusterEndcap().

double PFECALSuperClusterAlgo::threshPFClusterES_
private

Definition at line 142 of file PFECALSuperClusterAlgo.h.

Referenced by setThreshPFClusterES().

double PFECALSuperClusterAlgo::threshPFClusterSeed_
private

Definition at line 133 of file PFECALSuperClusterAlgo.h.

double PFECALSuperClusterAlgo::threshPFClusterSeedBarrel_
private

Definition at line 138 of file PFECALSuperClusterAlgo.h.

Referenced by setThreshPFClusterSeedBarrel().

double PFECALSuperClusterAlgo::threshPFClusterSeedEndcap_
private

Definition at line 140 of file PFECALSuperClusterAlgo.h.

Referenced by setThreshPFClusterSeedEndcap().

bool PFECALSuperClusterAlgo::usePS
private

Definition at line 156 of file PFECALSuperClusterAlgo.h.

Referenced by setUsePS().

bool PFECALSuperClusterAlgo::verbose_
private

Definition at line 131 of file PFECALSuperClusterAlgo.h.

Referenced by setVerbosityLevel().