CMS 3D CMS Logo

PFECALSuperClusterAlgo.cc
Go to the documentation of this file.
13 #include "Math/GenVector/VectorUtil.h"
14 #include "TFile.h"
15 #include "TH2F.h"
16 #include "TROOT.h"
17 #include "TMath.h"
18 
20 
21 #include <stdexcept>
22 #include <string>
23 #include <sstream>
24 #include <cmath>
25 
26 using namespace std;
27 namespace MK = reco::MustacheKernel;
28 
29 namespace {
30  typedef edm::View<reco::PFCluster> PFClusterView;
31  typedef edm::Ptr<reco::PFCluster> PFClusterPtr;
32  typedef edm::PtrVector<reco::PFCluster> PFClusterPtrVector;
33  typedef PFECALSuperClusterAlgo::CalibratedClusterPtr CalibClusterPtr;
34  typedef PFECALSuperClusterAlgo::CalibratedClusterPtrVector CalibClusterPtrVector;
35  typedef std::pair<reco::CaloClusterPtr::key_type,reco::CaloClusterPtr> EEPSPair;
36  typedef std::binary_function<const CalibClusterPtr&,
37  const CalibClusterPtr&,
38  bool> ClusBinaryFunction;
39 
40  typedef std::unary_function<const CalibClusterPtr&,
41  bool> ClusUnaryFunction;
42 
43  bool sortByKey(const EEPSPair& a, const EEPSPair& b) {
44  return a.first < b.first;
45  }
46 
47  inline double getPFClusterEnergy(const PFClusterPtr& p) {
48  return p->energy();
49  }
50 
51  inline double ptFast( const double energy,
52  const math::XYZPoint& position,
53  const math::XYZPoint& origin ) {
54  const auto v = position - origin;
55  return energy*std::sqrt(v.perp2()/v.mag2());
56  }
57 
58  struct SumPSEnergy : public std::binary_function<double,
59  const PFClusterPtr&,
60  double> {
61  PFLayer::Layer _thelayer;
62  SumPSEnergy(PFLayer::Layer layer) : _thelayer(layer) {}
63  double operator()(double a,
64  const PFClusterPtr& b) {
65  return a + (_thelayer == b->layer())*b->energy();
66  }
67  };
68 
69  struct GreaterByE : public ClusBinaryFunction {
70  bool operator()(const CalibClusterPtr& x,
71  const CalibClusterPtr& y) {
72  return x->energy() > y->energy() ;
73  }
74  };
75 
76  struct GreaterByEt : public ClusBinaryFunction {
77  bool operator()(const CalibClusterPtr& x,
78  const CalibClusterPtr& y) {
79  const math::XYZPoint zero(0,0,0);
80  const double xpt = ptFast(x->energy(),x->the_ptr()->position(),zero);
81  const double ypt = ptFast(y->energy(),y->the_ptr()->position(),zero);
82  return xpt > ypt;
83  }
84  };
85 
86  struct IsASeed : public ClusUnaryFunction {
87  const double threshold;
88  const bool cutET;
89  IsASeed(double thresh, bool useETcut = false) :
90  threshold(thresh), cutET(useETcut) {}
91  bool operator()(const CalibClusterPtr& x) {
92  const math::XYZPoint zero(0,0,0);
93  double e_or_et = x->energy();
94  if( cutET ) e_or_et = ptFast(e_or_et,x->the_ptr()->position(),zero);
95  return e_or_et > threshold;
96  }
97  };
98 
99  struct IsLinkedByRecHit : public ClusUnaryFunction {
100  const CalibClusterPtr the_seed;
101  const double _threshold, _majority;
102  const double _maxSatelliteDEta, _maxSatelliteDPhi;
103  double x_rechits_tot=0;
104  double x_rechits_match=0;
105  IsLinkedByRecHit(const CalibClusterPtr& s, const double threshold,
106  const double majority, const double maxDEta,
107  const double maxDPhi) :
108  the_seed(s),_threshold(threshold),_majority(majority),
109  _maxSatelliteDEta(maxDEta), _maxSatelliteDPhi(maxDPhi) {}
110  bool operator()(const CalibClusterPtr& x) {
111  if( the_seed->energy_nocalib() < _threshold ) return false;
112  const double dEta = std::abs(the_seed->eta()-x->eta());
113  const double dPhi =
114  std::abs(TVector2::Phi_mpi_pi(the_seed->phi() - x->phi()));
115  if( _maxSatelliteDEta < dEta || _maxSatelliteDPhi < dPhi) return false;
116  // now see if the clusters overlap in rechits
117  const auto& seedHitsAndFractions =
118  the_seed->the_ptr()->hitsAndFractions();
119  const auto& xHitsAndFractions =
120  x->the_ptr()->hitsAndFractions();
121  x_rechits_tot = xHitsAndFractions.size();
122  x_rechits_match = 0.0;
123  for( const std::pair<DetId, float>& seedHit : seedHitsAndFractions ) {
124  for( const std::pair<DetId, float>& xHit : xHitsAndFractions ) {
125  if( seedHit.first == xHit.first ) {
126  x_rechits_match += 1.0;
127  }
128  }
129  }
130  return x_rechits_match/x_rechits_tot > _majority;
131  }
132  };
133 
134  struct IsClustered : public ClusUnaryFunction {
135  const CalibClusterPtr the_seed;
137  bool dynamic_dphi;
138  double etawidthSuperCluster_ = .0 , phiwidthSuperCluster_ = .0;
139  IsClustered(const CalibClusterPtr s,
141  const bool dyn_dphi) :
142  the_seed(s), _type(ct), dynamic_dphi(dyn_dphi) {}
143  bool operator()(const CalibClusterPtr& x) {
144  const double dphi =
145  std::abs(TVector2::Phi_mpi_pi(the_seed->phi() - x->phi()));
146  const bool passes_dphi =
147  ( (!dynamic_dphi && dphi < phiwidthSuperCluster_ ) ||
148  (dynamic_dphi && MK::inDynamicDPhiWindow(the_seed->eta(),
149  the_seed->phi(),
150  x->energy_nocalib(),
151  x->eta(),
152  x->phi()) ) );
153 
154  switch( _type ) {
156  return ( std::abs(the_seed->eta()-x->eta())<etawidthSuperCluster_ &&
157  passes_dphi );
158  break;
160  return ( passes_dphi &&
161  MK::inMustache(the_seed->eta(),
162  the_seed->phi(),
163  x->energy_nocalib(),
164  x->eta(),
165  x->phi() ));
166  break;
167  default:
168  return false;
169  }
170  return false;
171  }
172  };
173 }
174 
176 
178 setPFClusterCalibration(const std::shared_ptr<PFEnergyCalibration>& calib) {
180 }
181 
183 
185  cc.consumes<edm::View<reco::PFCluster> >(iConfig.getParameter<edm::InputTag>("PFClusters"));
187  cc.consumes<reco::PFCluster::EEtoPSAssociation>(iConfig.getParameter<edm::InputTag>("ESAssociation"));
189  cc.consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>("BeamSpot"));
190 
191  if (useRegression_) {
192  const edm::ParameterSet &regconf = iConfig.getParameter<edm::ParameterSet>("regressionConfig");
193 
194  regr_.reset(new SCEnergyCorrectorSemiParm());
195  regr_->setTokens(regconf, cc);
196  }
197 
198  if (isOOTCollection_) { // OOT photons only
200  cc.consumes<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>("barrelRecHits"));
202  cc.consumes<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>("endcapRecHits"));
203  }
204 }
205 
207 
208  if (useRegression_) {
209  regr_->setEventSetup(setup);
210  }
211 
212  edm::ESHandle<ESEEIntercalibConstants> esEEInterCalibHandle_;
213  setup.get<ESEEIntercalibConstantsRcd>().get(esEEInterCalibHandle_);
214  _pfEnergyCalibration->initAlphaGamma_ESplanes_fromDB(esEEInterCalibHandle_.product());
215 
216  edm::ESHandle<ESChannelStatus> esChannelStatusHandle_;
217  setup.get<ESChannelStatusRcd>().get(esChannelStatusHandle_);
218  channelStatus_ = esChannelStatusHandle_.product();
219 
220 }
221 
222 
223 
226 
227  //load input collections
228  //Load the pfcluster collections
229  edm::Handle<edm::View<reco::PFCluster> > pfclustersHandle;
230  iEvent.getByToken( inputTagPFClusters_, pfclustersHandle );
231 
233  iEvent.getByToken( inputTagPFClustersES_, psAssociationHandle);
234 
235  const PFClusterView& clusters = *pfclustersHandle.product();
236  const reco::PFCluster::EEtoPSAssociation& psclusters = *psAssociationHandle.product();
237 
238  //load BeamSpot
240  iEvent.getByToken( inputTagBeamSpot_, bsHandle);
241  beamSpot_ = bsHandle.product();
242 
243  //initialize regression for this event
244  if (useRegression_) {
245  regr_->setEvent(iEvent);
246  }
247 
248  // reset the system for running
249  superClustersEB_ = std::make_unique<reco::SuperClusterCollection>();
250  _clustersEB.clear();
251  superClustersEE_ = std::make_unique<reco::SuperClusterCollection>();
252  _clustersEE.clear();
253  EEtoPS_ = &psclusters;
254 
255  //Select PF clusters available for the clustering
256  for ( size_t i = 0; i < clusters.size(); ++i ){
257  auto cluster = clusters.ptrAt(i);
258  LogDebug("PFClustering")
259  << "Loading PFCluster i="<<cluster.key()
260  <<" energy="<<cluster->energy()<<std::endl;
261 
262  // protection for sim clusters
263  if( cluster->caloID().detectors() == 0 &&
264  cluster->hitsAndFractions().empty() ) continue;
265 
266  CalibratedClusterPtr calib_cluster(new CalibratedPFCluster(cluster));
267  switch( cluster->layer() ) {
269  if( calib_cluster->energy() > threshPFClusterBarrel_ ) {
270  _clustersEB.push_back(calib_cluster);
271  }
272  break;
273  case PFLayer::HGCAL:
274  case PFLayer::ECAL_ENDCAP:
275  if( calib_cluster->energy() > threshPFClusterEndcap_ ) {
276  _clustersEE.push_back(calib_cluster);
277  }
278  break;
279  default:
280  break;
281  }
282  }
283  // sort full cluster collections by their calibrated energy
284  // this will put all the seeds first by construction
285  GreaterByEt greater;
286  std::sort(_clustersEB.begin(), _clustersEB.end(), greater);
287  std::sort(_clustersEE.begin(), _clustersEE.end(), greater);
288 
289  // set recHit collections for OOT photons
290  if (isOOTCollection_)
291  {
292  edm::Handle<EcalRecHitCollection> barrelRecHitsHandle;
293  iEvent.getByToken(inputTagBarrelRecHits_, barrelRecHitsHandle);
294  if (!barrelRecHitsHandle.isValid()) {
295  throw cms::Exception("PFECALSuperClusterAlgo")
296  << "If you use OOT photons, need to specify proper barrel rec hit collection";
297  }
298  barrelRecHits_ = barrelRecHitsHandle.product();
299 
300  edm::Handle<EcalRecHitCollection> endcapRecHitsHandle;
301  iEvent.getByToken(inputTagEndcapRecHits_, endcapRecHitsHandle);
302  if (!endcapRecHitsHandle.isValid()) {
303  throw cms::Exception("PFECALSuperClusterAlgo")
304  << "If you use OOT photons, need to specify proper endcap rec hit collection";
305  }
306  endcapRecHits_ = endcapRecHitsHandle.product();
307  }
308 }
309 
311  // clusterize the EB
313  // clusterize the EE
315 }
316 
318 buildAllSuperClusters(CalibClusterPtrVector& clusters,
319  double seedthresh) {
320  IsASeed seedable(seedthresh,threshIsET_);
321  // make sure only seeds appear at the front of the list of clusters
322  std::stable_partition(clusters.begin(),clusters.end(),seedable);
323  // in each iteration we are working on a list that is already sorted
324  // in the cluster energy and remains so through each iteration
325  // NB: since clusters is sorted in loadClusters any_of has O(1)
326  // timing until you run out of seeds!
327  while( std::any_of(clusters.cbegin(), clusters.cend(), seedable) ) {
328  buildSuperCluster(clusters.front(),clusters);
329  }
330 }
331 
333 buildSuperCluster(CalibClusterPtr& seed,
334  CalibClusterPtrVector& clusters) {
335  IsClustered IsClusteredWithSeed(seed,_clustype,_useDynamicDPhi);
336  IsLinkedByRecHit MatchesSeedByRecHit(seed,satelliteThreshold_,
337  fractionForMajority_,0.1,0.2);
338  bool isEE = false;
339  SumPSEnergy sumps1(PFLayer::PS1), sumps2(PFLayer::PS2);
340  switch( seed->the_ptr()->layer() ) {
342  IsClusteredWithSeed.phiwidthSuperCluster_ = phiwidthSuperClusterBarrel_;
343  IsClusteredWithSeed.etawidthSuperCluster_ = etawidthSuperClusterBarrel_;
344  edm::LogInfo("PFClustering") << "Building SC number "
345  << superClustersEB_->size() + 1
346  << " in the ECAL barrel!";
347  break;
348  case PFLayer::HGCAL:
349  case PFLayer::ECAL_ENDCAP:
350 
351  IsClusteredWithSeed.phiwidthSuperCluster_ = phiwidthSuperClusterEndcap_;
352  IsClusteredWithSeed.etawidthSuperCluster_ = etawidthSuperClusterEndcap_;
353  edm::LogInfo("PFClustering") << "Building SC number "
354  << superClustersEE_->size() + 1
355  << " in the ECAL endcap!" << std::endl;
356  isEE = true;
357  break;
358  default:
359  break;
360  }
361 
362  // this function shuffles the list of clusters into a list
363  // where all clustered sub-clusters are at the front
364  // and returns a pointer to the first unclustered cluster.
365  // The relative ordering of clusters is preserved
366  // (i.e. both resulting sub-lists are sorted by energy).
367  auto not_clustered = std::stable_partition(clusters.begin(),clusters.end(),
368  IsClusteredWithSeed);
369  // satellite cluster merging
370  // it was found that large clusters can split!
371  if( doSatelliteClusterMerge_ ) {
372  not_clustered = std::stable_partition(not_clustered,clusters.end(),
373  MatchesSeedByRecHit);
374  }
375 
376  if(verbose_) {
377  edm::LogInfo("PFClustering") << "Dumping cluster detail";
378  edm::LogVerbatim("PFClustering")
379  << "\tPassed seed: e = " << seed->energy_nocalib()
380  << " eta = " << seed->eta() << " phi = " << seed->phi()
381  << std::endl;
382  for( auto clus = clusters.cbegin(); clus != not_clustered; ++clus ) {
383  edm::LogVerbatim("PFClustering")
384  << "\t\tClustered cluster: e = " << (*clus)->energy_nocalib()
385  << " eta = " << (*clus)->eta() << " phi = " << (*clus)->phi()
386  << std::endl;
387  }
388  for( auto clus = not_clustered; clus != clusters.end(); ++clus ) {
389  edm::LogVerbatim("PFClustering")
390  << "\tNon-Clustered cluster: e = " << (*clus)->energy_nocalib()
391  << " eta = " << (*clus)->eta() << " phi = " << (*clus)->phi()
392  << std::endl;
393  }
394  }
395 
396  if (not_clustered == clusters.begin()) {
397  if(dropUnseedable_){
398  clusters.erase(clusters.begin());
399  return;
400  }
401  else {
402  throw cms::Exception("PFECALSuperClusterAlgo::buildSuperCluster")
403  << "Cluster is not seedable!" << std::endl
404  << "\tNon-Clustered cluster: e = " << (*not_clustered)->energy_nocalib()
405  << " eta = " << (*not_clustered)->eta() << " phi = " << (*not_clustered)->phi()
406  << std::endl;
407  }
408  }
409 
410  // move the clustered clusters out of available cluster list
411  // and into a temporary vector for building the SC
412  CalibratedClusterPtrVector clustered(clusters.begin(),not_clustered);
413  clusters.erase(clusters.begin(),not_clustered);
414  // need the vector of raw pointers for a PF width class
415  std::vector<const reco::PFCluster*> bare_ptrs;
416  // calculate necessary parameters and build the SC
417  double posX(0), posY(0), posZ(0),
418  corrSCEnergy(0), corrPS1Energy(0), corrPS2Energy(0),
419  ePS1(0), ePS2(0), energyweight(0), energyweighttot(0);
420  std::vector<double> ps1_energies, ps2_energies;
421  int condP1(1), condP2(1);
422  for( auto& clus : clustered ) {
423  ePS1 = ePS2 = 0;
424  energyweight = clus->energy_nocalib();
425  bare_ptrs.push_back(clus->the_ptr().get());
426  // update EE calibrated super cluster energies
427  if( isEE ) {
428  ePS1 = ePS2 = 0;
429  condP1 = condP2 = 1;
430  ps1_energies.clear();
431  ps2_energies.clear();
432  auto ee_key_val =
433  std::make_pair(clus->the_ptr().key(),edm::Ptr<reco::PFCluster>());
434  const auto clustops = std::equal_range(EEtoPS_->begin(),
435  EEtoPS_->end(),
436  ee_key_val,
437  sortByKey);
438  for( auto i_ps = clustops.first; i_ps != clustops.second; ++i_ps) {
439  edm::Ptr<reco::PFCluster> psclus(i_ps->second);
440 
441  auto const& recH_Frac = psclus->recHitFractions();
442 
443  switch( psclus->layer() ) {
444  case PFLayer::PS1:
445  ps1_energies.push_back(psclus->energy());
446  for (auto const& recH : recH_Frac){
447  ESDetId strip1 = recH.recHitRef()->detId();
448  if(strip1 != ESDetId(0)){
450  // getStatusCode() == 1 => dead channel
451  //apply correction if all recHits in dead region
452  if(status_p1->getStatusCode() == 0) condP1 = 0; //active
453  }
454  }
455  break;
456  case PFLayer::PS2:
457  ps2_energies.push_back(psclus->energy());
458  for (auto const& recH : recH_Frac){
459  ESDetId strip2 = recH.recHitRef()->detId();
460  if(strip2 != ESDetId(0)) {
462  if(status_p2->getStatusCode() == 0) condP2 = 0;
463  }
464  }
465  break;
466  default:
467  break;
468  }
469  }
470  if(condP1 == 1) ePS1 = -1.;
471  if(condP2 == 1) ePS2 = -1.;
472  _pfEnergyCalibration->energyEm(*(clus->the_ptr()),
473  ps1_energies,ps2_energies,
474  ePS1,ePS2,
476  }
477 
478  if(ePS1 == -1.) ePS1 = 0;
479  if(ePS2 == -1.) ePS2 = 0;
480 
481  switch( _eweight ) {
482  case kRaw: // energyweight is initialized to raw cluster energy
483  break;
484  case kCalibratedNoPS:
485  energyweight = clus->energy() - ePS1 - ePS2;
486  break;
487  case kCalibratedTotal:
488  energyweight = clus->energy();
489  break;
490  default:
491  break;
492  }
493  const math::XYZPoint& cluspos = clus->the_ptr()->position();
494  posX += energyweight * cluspos.X();
495  posY += energyweight * cluspos.Y();
496  posZ += energyweight * cluspos.Z();
497 
498  energyweighttot += energyweight;
499  corrSCEnergy += clus->energy();
500  corrPS1Energy += ePS1;
501  corrPS2Energy += ePS2;
502  }
503  posX /= energyweighttot;
504  posY /= energyweighttot;
505  posZ /= energyweighttot;
506 
507  // now build the supercluster
508  reco::SuperCluster new_sc(corrSCEnergy,math::XYZPoint(posX,posY,posZ));
509  new_sc.setCorrectedEnergy(corrSCEnergy);
510  new_sc.setSeed(clustered.front()->the_ptr());
511  new_sc.setPreshowerEnergy(corrPS1Energy+corrPS2Energy);
512  new_sc.setPreshowerEnergyPlane1(corrPS1Energy);
513  new_sc.setPreshowerEnergyPlane2(corrPS2Energy);
514  for( const auto& clus : clustered ) {
515  new_sc.addCluster(clus->the_ptr());
516 
517  auto& hits_and_fractions = clus->the_ptr()->hitsAndFractions();
518  for( auto& hit_and_fraction : hits_and_fractions ) {
519  new_sc.addHitAndFraction(hit_and_fraction.first,hit_and_fraction.second);
520  }
521  if( isEE ) {
522  auto ee_key_val =
523  std::make_pair(clus->the_ptr().key(),edm::Ptr<reco::PFCluster>());
524  const auto clustops = std::equal_range(EEtoPS_->begin(),
525  EEtoPS_->end(),
526  ee_key_val,
527  sortByKey);
528  // EE rechits should be uniquely matched to sets of pre-shower
529  // clusters at this point, so we throw an exception if otherwise
530  // now wrapped in EDM debug flags
531  for( auto i_ps = clustops.first; i_ps != clustops.second; ++i_ps) {
532  edm::Ptr<reco::PFCluster> psclus(i_ps->second);
533 #ifdef EDM_ML_DEBUG
534  auto found_pscluster = std::find_if(new_sc.preshowerClustersBegin(),
535  new_sc.preshowerClustersEnd(),
536  [&i_ps](const auto& i) { return i.key() == i_ps->first;});
537  if( found_pscluster == new_sc.preshowerClustersEnd() ) {
538 #endif
539  new_sc.addPreshowerCluster(psclus);
540 #ifdef EDM_ML_DEBUG
541  } else {
542  throw cms::Exception("PFECALSuperClusterAlgo::buildSuperCluster")
543  << "Found a PS cluster matched to more than one EE cluster!"
544  << std::endl << std::hex << psclus.get() << " == "
545  << found_pscluster->get() << std::dec << std::endl;
546  }
547 #endif
548  }
549  }
550  }
551 
552  // calculate linearly weighted cluster widths
553  PFClusterWidthAlgo pfwidth(bare_ptrs);
554  new_sc.setEtaWidth(pfwidth.pflowEtaWidth());
555  new_sc.setPhiWidth(pfwidth.pflowPhiWidth());
556 
557  // cache the value of the raw energy
558  new_sc.rawEnergy();
559 
560  //apply regression energy corrections
561  if( useRegression_ ) {
562  regr_->modifyObject(new_sc);
563  }
564 
565  // save the super cluster to the appropriate list (if it passes the final
566  // Et threshold)
567  //Note that Et is computed here with respect to the beamspot position
568  //in order to be consistent with the cut applied in the
569  //ElectronSeedProducer
570  double scEtBS =
571  ptFast(new_sc.energy(),new_sc.position(),beamSpot_->position());
572 
573  if ( scEtBS > threshSuperClusterEt_ ) {
574  switch( seed->the_ptr()->layer() ) {
576  if(isOOTCollection_) {
577  DetId seedId = new_sc.seed()->seed();
579  if (!seedRecHit->checkFlag(EcalRecHit::kOutOfTime)) break;
580  }
581  superClustersEB_->push_back(new_sc);
582  break;
583  case PFLayer::HGCAL:
584  case PFLayer::ECAL_ENDCAP:
585  if(isOOTCollection_) {
586  DetId seedId = new_sc.seed()->seed();
588  if (!seedRecHit->checkFlag(EcalRecHit::kOutOfTime)) break;
589  }
590  superClustersEE_->push_back(new_sc);
591  break;
592  default:
593  break;
594  }
595  }
596 }
#define LogDebug(id)
PFLayer::Layer layer() const
cluster layer, see PFLayer.h in this directory
Definition: PFCluster.cc:128
T getParameter(std::string const &) const
const math::XYZPoint & position() const
cluster centroid position
Definition: CaloCluster.h:131
CaloCluster_iterator preshowerClustersBegin() const
fist iterator over PreshowerCluster constituents
Definition: SuperCluster.h:81
const ESChannelStatus * channelStatus_
void addHitAndFraction(DetId id, float fraction)
Definition: CaloCluster.h:190
edm::EDGetTokenT< edm::View< reco::PFCluster > > inputTagPFClusters_
std::unique_ptr< SCEnergyCorrectorSemiParm > regr_
bool inDynamicDPhiWindow(const float seedEta, const float seedPhi, const float ClustE, const float ClusEta, const float clusPhi)
Definition: Mustache.cc:73
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
std::unique_ptr< reco::SuperClusterCollection > superClustersEB_
void setPreshowerEnergyPlane2(double preshowerEnergy2)
Definition: SuperCluster.h:61
edm::EDGetTokenT< EcalRecHitCollection > inputTagBarrelRecHits_
const self & getMap() const
T const * get() const
Returns C++ pointer to the item.
Definition: Ptr.h:159
double pflowPhiWidth() const
std::vector< EcalRecHit >::const_iterator const_iterator
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:1
const reco::BeamSpot * beamSpot_
const reco::PFCluster::EEtoPSAssociation * EEtoPS_
void setSeed(const CaloClusterPtr &r)
list of used xtals by DetId // now inherited by CaloCluster
Definition: SuperCluster.h:96
#define nullptr
void setPhiWidth(double pw)
Definition: SuperCluster.h:62
double pflowEtaWidth() const
double ptFast(const double energy, const math::XYZPoint &position, const math::XYZPoint &origin)
std::unique_ptr< reco::SuperClusterCollection > superClustersEE_
edm::EDGetTokenT< EcalRecHitCollection > inputTagEndcapRecHits_
void loadAndSortPFClusters(const edm::Event &evt)
void update(const edm::EventSetup &)
void setEtaWidth(double ew)
Definition: SuperCluster.h:63
std::vector< CalibratedClusterPtr > CalibratedClusterPtrVector
std::shared_ptr< CalibratedPFCluster > CalibratedClusterPtr
MVATrainerComputer * calib
Definition: MVATrainer.cc:64
int iEvent
Definition: GenABIO.cc:230
void buildAllSuperClusters(CalibratedClusterPtrVector &, double seedthresh)
std::shared_ptr< PFEnergyCalibration > _pfEnergyCalibration
void setTokens(const edm::ParameterSet &, edm::ConsumesCollector &&)
void setCorrectedEnergy(double cenergy)
Definition: CaloCluster.h:114
T sqrt(T t)
Definition: SSEVec.h:18
const_iterator find(uint32_t rawId) const
std::vector< std::pair< CaloClusterPtr::key_type, edm::Ptr< PFCluster > > > EEtoPSAssociation
Definition: PFCluster.h:50
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
CalibratedClusterPtrVector _clustersEE
double energy() const
cluster energy
Definition: CaloCluster.h:126
bool isValid() const
Definition: HandleBase.h:74
double energy() const
cluster energy
Definition: PFCluster.h:82
Layer
layer definition
Definition: PFLayer.h:31
Definition: DetId.h:18
double rawEnergy() const
raw uncorrected energy (sum of energies of component BasicClusters)
Definition: SuperCluster.h:47
edm::EDGetTokenT< reco::PFCluster::EEtoPSAssociation > inputTagPFClustersES_
T const * product() const
Definition: Handle.h:81
edm::EDGetTokenT< reco::BeamSpot > inputTagBeamSpot_
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
reco::PFCluster::EEtoPSAssociation::value_type EEPSPair
void addPreshowerCluster(const CaloClusterPtr &r)
add reference to constituent BasicCluster
Definition: SuperCluster.h:117
bool sortByKey(const EEPSPair &a, const EEPSPair &b)
double b
Definition: hdecay.h:120
void buildSuperCluster(CalibratedClusterPtr &, CalibratedClusterPtrVector &)
std::vector< Item >::const_iterator const_iterator
bool GreaterByE(const T &a1, const T &a2)
void addCluster(const CaloClusterPtr &r)
add reference to constituent BasicCluster
Definition: SuperCluster.h:111
iterator find(key_type k)
double a
Definition: hdecay.h:121
static int position[264][3]
Definition: ReadPGInfo.cc:509
T get() const
Definition: EventSetup.h:63
const std::vector< reco::PFRecHitFraction > & recHitFractions() const
vector of rechit fractions
Definition: PFCluster.h:72
const Point & position() const
position
Definition: BeamSpot.h:62
const CaloClusterPtr & seed() const
seed BasicCluster
Definition: SuperCluster.h:66
const EcalRecHitCollection * barrelRecHits_
T const * product() const
Definition: ESHandle.h:86
void setPreshowerEnergyPlane1(double preshowerEnergy1)
Definition: SuperCluster.h:60
void setPFClusterCalibration(const std::shared_ptr< PFEnergyCalibration > &)
CalibratedClusterPtrVector _clustersEB
const EcalRecHitCollection * endcapRecHits_
bool inMustache(const float maxEta, const float maxPhi, const float ClustE, const float ClusEta, const float ClusPhi)
Definition: Mustache.cc:9
CaloCluster_iterator preshowerClustersEnd() const
last iterator over PreshowerCluster constituents
Definition: SuperCluster.h:84
void setPreshowerEnergy(double preshowerEnergy)
Definition: SuperCluster.h:59