CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFECALSuperClusterAlgo.cc
Go to the documentation of this file.
7 #include "Math/GenVector/VectorUtil.h"
8 #include "TFile.h"
9 #include "TH2F.h"
10 #include "TROOT.h"
11 #include "TMath.h"
12 
14 
15 #include <stdexcept>
16 #include <string>
17 #include <sstream>
18 #include <cmath>
19 
20 using namespace std;
21 namespace MK = reco::MustacheKernel;
22 
23 namespace {
24  typedef edm::View<reco::PFCluster> PFClusterView;
25  typedef edm::Ptr<reco::PFCluster> PFClusterPtr;
26  typedef edm::PtrVector<reco::PFCluster> PFClusterPtrVector;
27  typedef PFECALSuperClusterAlgo::CalibratedClusterPtr CalibClusterPtr;
28  typedef PFECALSuperClusterAlgo::CalibratedClusterPtrVector CalibClusterPtrVector;
29  typedef std::pair<reco::CaloClusterPtr::key_type,reco::CaloClusterPtr> EEPSPair;
30  typedef std::binary_function<const CalibClusterPtr&,
31  const CalibClusterPtr&,
32  bool> ClusBinaryFunction;
33 
34  typedef std::unary_function<const CalibClusterPtr&,
35  bool> ClusUnaryFunction;
36 
37  bool sortByKey(const EEPSPair& a, const EEPSPair& b) {
38  return a.first < b.first;
39  }
40 
41  inline double getPFClusterEnergy(const PFClusterPtr& p) {
42  return p->energy();
43  }
44 
45  inline double ptFast( const double energy,
46  const math::XYZPoint& position,
47  const math::XYZPoint& origin ) {
48  const auto v = position - origin;
49  return energy*std::sqrt(v.perp2()/v.mag2());
50  }
51 
52  struct SumPSEnergy : public std::binary_function<double,
53  const PFClusterPtr&,
54  double> {
55  PFLayer::Layer _thelayer;
56  SumPSEnergy(PFLayer::Layer layer) : _thelayer(layer) {}
57  double operator()(double a,
58  const PFClusterPtr& b) {
59  return a + (_thelayer == b->layer())*b->energy();
60  }
61  };
62 
63  struct GreaterByE : public ClusBinaryFunction {
64  bool operator()(const CalibClusterPtr& x,
65  const CalibClusterPtr& y) {
66  return x->energy() > y->energy() ;
67  }
68  };
69 
70  struct GreaterByEt : public ClusBinaryFunction {
71  bool operator()(const CalibClusterPtr& x,
72  const CalibClusterPtr& y) {
73  const math::XYZPoint zero(0,0,0);
74  const double xpt = ptFast(x->energy(),x->the_ptr()->position(),zero);
75  const double ypt = ptFast(y->energy(),y->the_ptr()->position(),zero);
76  return xpt > ypt;
77  }
78  };
79 
80  struct IsASeed : public ClusUnaryFunction {
81  const double threshold;
82  const bool cutET;
83  IsASeed(double thresh, bool useETcut = false) :
84  threshold(thresh), cutET(useETcut) {}
85  bool operator()(const CalibClusterPtr& x) {
86  const math::XYZPoint zero(0,0,0);
87  double e_or_et = x->energy();
88  if( cutET ) e_or_et = ptFast(e_or_et,x->the_ptr()->position(),zero);
89  return e_or_et > threshold;
90  }
91  };
92 
93  struct IsLinkedByRecHit : public ClusUnaryFunction {
94  const CalibClusterPtr the_seed;
95  const double _threshold, _majority;
96  const double _maxSatelliteDEta, _maxSatelliteDPhi;
97  double x_rechits_tot=0;
98  double x_rechits_match=0;
99  IsLinkedByRecHit(const CalibClusterPtr& s, const double threshold,
100  const double majority, const double maxDEta,
101  const double maxDPhi) :
102  the_seed(s),_threshold(threshold),_majority(majority),
103  _maxSatelliteDEta(maxDEta), _maxSatelliteDPhi(maxDPhi) {}
104  bool operator()(const CalibClusterPtr& x) {
105  if( the_seed->energy_nocalib() < _threshold ) return false;
106  const double dEta = std::abs(the_seed->eta()-x->eta());
107  const double dPhi =
108  std::abs(TVector2::Phi_mpi_pi(the_seed->phi() - x->phi()));
109  if( _maxSatelliteDEta < dEta || _maxSatelliteDPhi < dPhi) return false;
110  // now see if the clusters overlap in rechits
111  const auto& seedHitsAndFractions =
112  the_seed->the_ptr()->hitsAndFractions();
113  const auto& xHitsAndFractions =
114  x->the_ptr()->hitsAndFractions();
115  x_rechits_tot = xHitsAndFractions.size();
116  x_rechits_match = 0.0;
117  for( const std::pair<DetId, float>& seedHit : seedHitsAndFractions ) {
118  for( const std::pair<DetId, float>& xHit : xHitsAndFractions ) {
119  if( seedHit.first == xHit.first ) {
120  x_rechits_match += 1.0;
121  }
122  }
123  }
124  return x_rechits_match/x_rechits_tot > _majority;
125  }
126  };
127 
128  struct IsClustered : public ClusUnaryFunction {
129  const CalibClusterPtr the_seed;
131  bool dynamic_dphi;
132  double etawidthSuperCluster_ = .0 , phiwidthSuperCluster_ = .0;
133  IsClustered(const CalibClusterPtr s,
135  const bool dyn_dphi) :
136  the_seed(s), _type(ct), dynamic_dphi(dyn_dphi) {}
137  bool operator()(const CalibClusterPtr& x) {
138  const double dphi =
139  std::abs(TVector2::Phi_mpi_pi(the_seed->phi() - x->phi()));
140  const bool passes_dphi =
141  ( (!dynamic_dphi && dphi < phiwidthSuperCluster_ ) ||
142  (dynamic_dphi && MK::inDynamicDPhiWindow(the_seed->eta(),
143  the_seed->phi(),
144  x->energy_nocalib(),
145  x->eta(),
146  x->phi()) ) );
147 
148  switch( _type ) {
150  return ( std::abs(the_seed->eta()-x->eta())<etawidthSuperCluster_ &&
151  passes_dphi );
152  break;
154  return ( passes_dphi &&
155  MK::inMustache(the_seed->eta(),
156  the_seed->phi(),
157  x->energy_nocalib(),
158  x->eta(),
159  x->phi() ));
160  break;
161  default:
162  return false;
163  }
164  return false;
165  }
166  };
167 }
168 
170 
172 setPFClusterCalibration(const std::shared_ptr<PFEnergyCalibration>& calib) {
174 }
175 
177 
179  cc.consumes<edm::View<reco::PFCluster> >(iConfig.getParameter<edm::InputTag>("PFClusters"));
181  cc.consumes<reco::PFCluster::EEtoPSAssociation>(iConfig.getParameter<edm::InputTag>("ESAssociation"));
183  cc.consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>("BeamSpot"));
184 
185  if (useRegression_) {
186  const edm::ParameterSet &regconf = iConfig.getParameter<edm::ParameterSet>("regressionConfig");
187 
188  regr_.reset(new PFSCRegressionCalc(regconf));
189  regr_->varCalc()->setTokens(regconf,std::move(cc));
190  }
191 
192 }
193 
195 
196  if (useRegression_) {
197  regr_->update(setup);
198  }
199 
200 }
201 
202 
203 
206 
207  //load input collections
208  //Load the pfcluster collections
209  edm::Handle<edm::View<reco::PFCluster> > pfclustersHandle;
210  iEvent.getByToken( inputTagPFClusters_, pfclustersHandle );
211 
213  iEvent.getByToken( inputTagPFClustersES_, psAssociationHandle);
214 
215  const PFClusterView& clusters = *pfclustersHandle.product();
216  const reco::PFCluster::EEtoPSAssociation& psclusters = *psAssociationHandle.product();
217 
218  //load BeamSpot
220  iEvent.getByToken( inputTagBeamSpot_, bsHandle);
221  beamSpot_ = bsHandle.product();
222 
223  //initialize regression for this event
224  if (useRegression_) {
225  regr_->varCalc()->setEvent(iEvent);
226  }
227 
228  // reset the system for running
230  _clustersEB.clear();
232  _clustersEE.clear();
233  EEtoPS_ = &psclusters;
234 
235  //Select PF clusters available for the clustering
236  for ( size_t i = 0; i < clusters.size(); ++i ){
237  auto cluster = clusters.ptrAt(i);
238  LogDebug("PFClustering")
239  << "Loading PFCluster i="<<cluster.key()
240  <<" energy="<<cluster->energy()<<std::endl;
241 
242  CalibratedClusterPtr calib_cluster(new CalibratedPFCluster(cluster));
243  switch( cluster->layer() ) {
245  if( calib_cluster->energy() > threshPFClusterBarrel_ ) {
246  _clustersEB.push_back(calib_cluster);
247  }
248  break;
250  if( calib_cluster->energy() > threshPFClusterEndcap_ ) {
251  _clustersEE.push_back(calib_cluster);
252  }
253  break;
254  default:
255  break;
256  }
257  }
258  // sort full cluster collections by their calibrated energy
259  // this will put all the seeds first by construction
260  GreaterByEt greater;
261  std::sort(_clustersEB.begin(), _clustersEB.end(), greater);
262  std::sort(_clustersEE.begin(), _clustersEE.end(), greater);
263 
264 }
265 
267  // clusterize the EB
269  // clusterize the EE
271 }
272 
274 buildAllSuperClusters(CalibClusterPtrVector& clusters,
275  double seedthresh) {
276  IsASeed seedable(seedthresh,threshIsET_);
277  // make sure only seeds appear at the front of the list of clusters
278  std::stable_partition(clusters.begin(),clusters.end(),seedable);
279  // in each iteration we are working on a list that is already sorted
280  // in the cluster energy and remains so through each iteration
281  // NB: since clusters is sorted in loadClusters any_of has O(1)
282  // timing until you run out of seeds!
283  while( std::any_of(clusters.cbegin(), clusters.cend(), seedable) ) {
284  buildSuperCluster(clusters.front(),clusters);
285  }
286 }
287 
289 buildSuperCluster(CalibClusterPtr& seed,
290  CalibClusterPtrVector& clusters) {
291  IsClustered IsClusteredWithSeed(seed,_clustype,_useDynamicDPhi);
292  IsLinkedByRecHit MatchesSeedByRecHit(seed,satelliteThreshold_,
293  fractionForMajority_,0.1,0.2);
294  bool isEE = false;
295  SumPSEnergy sumps1(PFLayer::PS1), sumps2(PFLayer::PS2);
296  switch( seed->the_ptr()->layer() ) {
298  IsClusteredWithSeed.phiwidthSuperCluster_ = phiwidthSuperClusterBarrel_;
299  IsClusteredWithSeed.etawidthSuperCluster_ = etawidthSuperClusterBarrel_;
300  edm::LogInfo("PFClustering") << "Building SC number "
301  << superClustersEB_->size() + 1
302  << " in the ECAL barrel!";
303  break;
304  case PFLayer::ECAL_ENDCAP:
305  IsClusteredWithSeed.phiwidthSuperCluster_ = phiwidthSuperClusterEndcap_;
306  IsClusteredWithSeed.etawidthSuperCluster_ = etawidthSuperClusterEndcap_;
307  edm::LogInfo("PFClustering") << "Building SC number "
308  << superClustersEE_->size() + 1
309  << " in the ECAL endcap!" << std::endl;
310  isEE = true;
311  break;
312  default:
313  break;
314  }
315 
316  // this function shuffles the list of clusters into a list
317  // where all clustered sub-clusters are at the front
318  // and returns a pointer to the first unclustered cluster.
319  // The relative ordering of clusters is preserved
320  // (i.e. both resulting sub-lists are sorted by energy).
321  auto not_clustered = std::stable_partition(clusters.begin(),clusters.end(),
322  IsClusteredWithSeed);
323  // satellite cluster merging
324  // it was found that large clusters can split!
325  if( doSatelliteClusterMerge_ ) {
326  not_clustered = std::stable_partition(not_clustered,clusters.end(),
327  MatchesSeedByRecHit);
328  }
329 
330  if(verbose_) {
331  edm::LogInfo("PFClustering") << "Dumping cluster detail";
332  edm::LogVerbatim("PFClustering")
333  << "\tPassed seed: e = " << seed->energy_nocalib()
334  << " eta = " << seed->eta() << " phi = " << seed->phi()
335  << std::endl;
336  for( auto clus = clusters.cbegin(); clus != not_clustered; ++clus ) {
337  edm::LogVerbatim("PFClustering")
338  << "\t\tClustered cluster: e = " << (*clus)->energy_nocalib()
339  << " eta = " << (*clus)->eta() << " phi = " << (*clus)->phi()
340  << std::endl;
341  }
342  for( auto clus = not_clustered; clus != clusters.end(); ++clus ) {
343  edm::LogVerbatim("PFClustering")
344  << "\tNon-Clustered cluster: e = " << (*clus)->energy_nocalib()
345  << " eta = " << (*clus)->eta() << " phi = " << (*clus)->phi()
346  << std::endl;
347  }
348  }
349  // move the clustered clusters out of available cluster list
350  // and into a temporary vector for building the SC
351  CalibratedClusterPtrVector clustered(clusters.begin(),not_clustered);
352  clusters.erase(clusters.begin(),not_clustered);
353  // need the vector of raw pointers for a PF width class
354  std::vector<const reco::PFCluster*> bare_ptrs;
355  // calculate necessary parameters and build the SC
356  double posX(0), posY(0), posZ(0),
357  corrSCEnergy(0), corrPS1Energy(0), corrPS2Energy(0),
358  ePS1(0), ePS2(0), energyweight(0), energyweighttot(0);
359  std::vector<double> ps1_energies, ps2_energies;
360  for( auto& clus : clustered ) {
361  ePS1 = ePS2 = 0;
362  energyweight = clus->energy_nocalib();
363  bare_ptrs.push_back(clus->the_ptr().get());
364  // update EE calibrated super cluster energies
365  if( isEE ) {
366  ps1_energies.clear();
367  ps2_energies.clear();
368  auto ee_key_val =
369  std::make_pair(clus->the_ptr().key(),edm::Ptr<reco::PFCluster>());
370  const auto clustops = std::equal_range(EEtoPS_->begin(),
371  EEtoPS_->end(),
372  ee_key_val,
373  sortByKey);
374  for( auto i_ps = clustops.first; i_ps != clustops.second; ++i_ps) {
375  edm::Ptr<reco::PFCluster> psclus(i_ps->second);
376  switch( psclus->layer() ) {
377  case PFLayer::PS1:
378  ps1_energies.push_back(psclus->energy());
379  break;
380  case PFLayer::PS2:
381  ps2_energies.push_back(psclus->energy());
382  break;
383  default:
384  break;
385  }
386  }
387  _pfEnergyCalibration->energyEm(*(clus->the_ptr()),
388  ps1_energies,ps2_energies,
389  ePS1,ePS2,
391  }
392 
393  switch( _eweight ) {
394  case kRaw: // energyweight is initialized to raw cluster energy
395  break;
396  case kCalibratedNoPS:
397  energyweight = clus->energy() - ePS1 - ePS2;
398  break;
399  case kCalibratedTotal:
400  energyweight = clus->energy();
401  break;
402  default:
403  break;
404  }
405  const math::XYZPoint& cluspos = clus->the_ptr()->position();
406  posX += energyweight * cluspos.X();
407  posY += energyweight * cluspos.Y();
408  posZ += energyweight * cluspos.Z();
409 
410  energyweighttot += energyweight;
411  corrSCEnergy += clus->energy();
412  corrPS1Energy += ePS1;
413  corrPS2Energy += ePS2;
414  }
415  posX /= energyweighttot;
416  posY /= energyweighttot;
417  posZ /= energyweighttot;
418 
419  // now build the supercluster
420  reco::SuperCluster new_sc(corrSCEnergy,math::XYZPoint(posX,posY,posZ));
421  new_sc.setCorrectedEnergy(corrSCEnergy);
422  new_sc.setSeed(clustered.front()->the_ptr());
423  new_sc.setPreshowerEnergy(corrPS1Energy+corrPS2Energy);
424  new_sc.setPreshowerEnergyPlane1(corrPS1Energy);
425  new_sc.setPreshowerEnergyPlane2(corrPS2Energy);
426  for( const auto& clus : clustered ) {
427  new_sc.addCluster(clus->the_ptr());
428 
429  auto& hits_and_fractions = clus->the_ptr()->hitsAndFractions();
430  for( auto& hit_and_fraction : hits_and_fractions ) {
431  new_sc.addHitAndFraction(hit_and_fraction.first,hit_and_fraction.second);
432  }
433  if( isEE ) {
434  auto ee_key_val =
435  std::make_pair(clus->the_ptr().key(),edm::Ptr<reco::PFCluster>());
436  const auto clustops = std::equal_range(EEtoPS_->begin(),
437  EEtoPS_->end(),
438  ee_key_val,
439  sortByKey);
440  // EE rechits should be uniquely matched to sets of pre-shower
441  // clusters at this point, so we throw an exception if otherwise
442  // now wrapped in EDM debug flags
443  for( auto i_ps = clustops.first; i_ps != clustops.second; ++i_ps) {
444  edm::Ptr<reco::PFCluster> psclus(i_ps->second);
445 #ifdef EDM_ML_DEBUG
446  auto found_pscluster = std::find(new_sc.preshowerClustersBegin(),
447  new_sc.preshowerClustersEnd(),
448  i_ps->second);
449  if( found_pscluster == new_sc.preshowerClustersEnd() ) {
450 #endif
451  new_sc.addPreshowerCluster(psclus);
452 #ifdef EDM_ML_DEBUG
453  } else {
454  throw cms::Exception("PFECALSuperClusterAlgo::buildSuperCluster")
455  << "Found a PS cluster matched to more than one EE cluster!"
456  << std::endl << std::hex << psclus.get() << " == "
457  << found_pscluster->get() << std::dec << std::endl;
458  }
459 #endif
460  }
461  }
462  }
463 
464  // calculate linearly weighted cluster widths
465  PFClusterWidthAlgo pfwidth(bare_ptrs);
466  new_sc.setEtaWidth(pfwidth.pflowEtaWidth());
467  new_sc.setPhiWidth(pfwidth.pflowPhiWidth());
468 
469  // cache the value of the raw energy
470  new_sc.rawEnergy();
471 
472  //apply regression energy corrections
473  if( useRegression_ ) {
474  double cor = regr_->getCorrection(new_sc);
475  new_sc.setEnergy(cor*new_sc.correctedEnergy());
476  }
477 
478  // save the super cluster to the appropriate list (if it passes the final
479  // Et threshold)
480  //Note that Et is computed here with respect to the beamspot position
481  //in order to be consistent with the cut applied in the
482  //ElectronSeedProducer
483  double scEtBS =
484  ptFast(new_sc.energy(),new_sc.position(),beamSpot_->position());
485 
486  if ( scEtBS > threshSuperClusterEt_ ) {
487  switch( seed->the_ptr()->layer() ) {
489  superClustersEB_->push_back(new_sc);
490  break;
491  case PFLayer::ECAL_ENDCAP:
492  superClustersEE_->push_back(new_sc);
493  break;
494  default:
495  break;
496  }
497  }
498 }
#define LogDebug(id)
T getParameter(std::string const &) const
const math::XYZPoint & position() const
cluster centroid position
Definition: CaloCluster.h:126
CaloCluster_iterator preshowerClustersBegin() const
fist iterator over PreshowerCluster constituents
Definition: SuperCluster.h:81
std::auto_ptr< reco::SuperClusterCollection > superClustersEB_
int i
Definition: DBlmapReader.cc:9
double correctedEnergy() const
Definition: CaloCluster.h:122
void addHitAndFraction(DetId id, float fraction)
Definition: CaloCluster.h:185
edm::EDGetTokenT< edm::View< reco::PFCluster > > inputTagPFClusters_
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:464
void setPreshowerEnergyPlane2(double preshowerEnergy2)
Definition: SuperCluster.h:61
T const * get() const
Returns C++ pointer to the item.
Definition: Ptr.h:160
double pflowPhiWidth() const
const reco::BeamSpot * beamSpot_
void setEnergy(double energy)
Definition: CaloCluster.h:108
const reco::PFCluster::EEtoPSAssociation * EEtoPS_
double Phi_mpi_pi(double x)
Definition: JetUtil.h:24
void setSeed(const CaloClusterPtr &r)
list of used xtals by DetId // now inherited by CaloCluster
Definition: SuperCluster.h:96
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
void setPhiWidth(double pw)
Definition: SuperCluster.h:62
double pflowEtaWidth() const
double ptFast(const double energy, const math::XYZPoint &position, const math::XYZPoint &origin)
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 &&)
double dPhi(double phi1, double phi2)
Definition: JetUtil.h:30
std::vector< SuperCluster > SuperClusterCollection
collection of SuperCluser objectr
void setCorrectedEnergy(double cenergy)
Definition: CaloCluster.h:109
T sqrt(T t)
Definition: SSEVec.h:48
std::unique_ptr< PFSCRegressionCalc > regr_
def move
Definition: eostools.py:510
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:121
Layer
layer definition
Definition: PFLayer.h:31
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
void addPreshowerCluster(const CaloClusterPtr &r)
add reference to constituent BasicCluster
Definition: SuperCluster.h:117
std::auto_ptr< reco::SuperClusterCollection > superClustersEE_
double b
Definition: hdecay.h:120
void buildSuperCluster(CalibratedClusterPtr &, CalibratedClusterPtrVector &)
bool GreaterByE(const T &a1, const T &a2)
void addCluster(const CaloClusterPtr &r)
add reference to constituent BasicCluster
Definition: SuperCluster.h:111
double a
Definition: hdecay.h:121
static int position[264][3]
Definition: ReadPGInfo.cc:509
const Point & position() const
position
Definition: BeamSpot.h:62
void setPreshowerEnergyPlane1(double preshowerEnergy1)
Definition: SuperCluster.h:60
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")
void setPFClusterCalibration(const std::shared_ptr< PFEnergyCalibration > &)
CalibratedClusterPtrVector _clustersEB
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
SCRegressionCalculator< BaselinePFSCRegression > PFSCRegressionCalc
void setPreshowerEnergy(double preshowerEnergy)
Definition: SuperCluster.h:59