7 #include "Math/GenVector/VectorUtil.h"
21 namespace MK = reco::MustacheKernel;
29 typedef std::pair<reco::CaloClusterPtr::key_type,reco::CaloClusterPtr> EEPSPair;
30 typedef std::binary_function<
const CalibClusterPtr&,
31 const CalibClusterPtr&,
32 bool> ClusBinaryFunction;
34 typedef std::unary_function<
const CalibClusterPtr&,
35 bool> ClusUnaryFunction;
37 bool sortByKey(
const EEPSPair&
a,
const EEPSPair&
b) {
38 return a.first < b.first;
41 inline double getPFClusterEnergy(
const PFClusterPtr&
p) {
48 const auto v = position - origin;
52 struct SumPSEnergy :
public std::binary_function<double,
57 double operator()(
double a,
58 const PFClusterPtr&
b) {
59 return a + (_thelayer == b->layer())*b->energy();
63 struct GreaterByE :
public ClusBinaryFunction {
64 bool operator()(
const CalibClusterPtr&
x,
65 const CalibClusterPtr&
y) {
66 return x->energy() > y->energy() ;
70 struct GreaterByEt :
public ClusBinaryFunction {
71 bool operator()(
const CalibClusterPtr& x,
72 const CalibClusterPtr& y) {
74 const double xpt =
ptFast(x->energy(),x->the_ptr()->position(),zero);
75 const double ypt =
ptFast(y->energy(),y->the_ptr()->position(),zero);
80 struct IsASeed :
public ClusUnaryFunction {
83 IsASeed(
double thresh,
bool useETcut =
false) :
85 bool operator()(
const CalibClusterPtr& x) {
87 double e_or_et = x->energy();
88 if( cutET ) e_or_et =
ptFast(e_or_et,x->the_ptr()->position(),zero);
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());
109 if( _maxSatelliteDEta < dEta || _maxSatelliteDPhi < dPhi)
return false;
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;
124 return x_rechits_match/x_rechits_tot > _majority;
128 struct IsClustered :
public ClusUnaryFunction {
129 const CalibClusterPtr the_seed;
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) {
140 const bool passes_dphi =
141 ( (!dynamic_dphi && dphi < phiwidthSuperCluster_ ) ||
150 return (
std::abs(the_seed->eta()-x->eta())<etawidthSuperCluster_ &&
154 return ( passes_dphi &&
189 regr_->varCalc()->setTokens(regconf,std::move(cc));
197 regr_->update(setup);
215 const PFClusterView& clusters = *pfclustersHandle.
product();
225 regr_->varCalc()->setEvent(iEvent);
236 for (
size_t i = 0;
i < clusters.size(); ++
i ){
237 auto cluster = clusters.ptrAt(
i);
239 <<
"Loading PFCluster i="<<cluster.key()
240 <<
" energy="<<cluster->energy()<<std::endl;
243 switch( cluster->layer() ) {
278 std::stable_partition(clusters.begin(),clusters.end(),seedable);
283 while( std::any_of(clusters.cbegin(), clusters.cend(), seedable) ) {
290 CalibClusterPtrVector& clusters) {
296 switch( seed->the_ptr()->layer() ) {
302 <<
" in the ECAL barrel!";
309 <<
" in the ECAL endcap!" << std::endl;
321 auto not_clustered = std::stable_partition(clusters.begin(),clusters.end(),
322 IsClusteredWithSeed);
326 not_clustered = std::stable_partition(not_clustered,clusters.end(),
327 MatchesSeedByRecHit);
331 edm::LogInfo(
"PFClustering") <<
"Dumping cluster detail";
333 <<
"\tPassed seed: e = " << seed->energy_nocalib()
334 <<
" eta = " << seed->eta() <<
" phi = " << seed->phi()
336 for(
auto clus = clusters.cbegin(); clus != not_clustered; ++clus ) {
338 <<
"\t\tClustered cluster: e = " << (*clus)->energy_nocalib()
339 <<
" eta = " << (*clus)->eta() <<
" phi = " << (*clus)->phi()
342 for(
auto clus = not_clustered; clus != clusters.end(); ++clus ) {
344 <<
"\tNon-Clustered cluster: e = " << (*clus)->energy_nocalib()
345 <<
" eta = " << (*clus)->eta() <<
" phi = " << (*clus)->phi()
352 clusters.erase(clusters.begin(),not_clustered);
354 std::vector<const reco::PFCluster*> bare_ptrs;
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 ) {
362 energyweight = clus->energy_nocalib();
363 bare_ptrs.push_back(clus->the_ptr().get());
366 ps1_energies.clear();
367 ps2_energies.clear();
370 const auto clustops = std::equal_range(
EEtoPS_->begin(),
374 for(
auto i_ps = clustops.first; i_ps != clustops.second; ++i_ps) {
376 switch( psclus->layer() ) {
378 ps1_energies.push_back(psclus->energy());
381 ps2_energies.push_back(psclus->energy());
388 ps1_energies,ps2_energies,
397 energyweight = clus->energy() - ePS1 - ePS2;
400 energyweight = clus->energy();
406 posX += energyweight * cluspos.X();
407 posY += energyweight * cluspos.Y();
408 posZ += energyweight * cluspos.Z();
410 energyweighttot += energyweight;
411 corrSCEnergy += clus->energy();
412 corrPS1Energy += ePS1;
413 corrPS2Energy += ePS2;
415 posX /= energyweighttot;
416 posY /= energyweighttot;
417 posZ /= energyweighttot;
422 new_sc.
setSeed(clustered.front()->the_ptr());
426 for(
const auto& clus : clustered ) {
429 auto& hits_and_fractions = clus->the_ptr()->hitsAndFractions();
430 for(
auto& hit_and_fraction : hits_and_fractions ) {
436 const auto clustops = std::equal_range(
EEtoPS_->begin(),
443 for(
auto i_ps = clustops.first; i_ps != clustops.second; ++i_ps) {
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;
474 double cor =
regr_->getCorrection(new_sc);
487 switch( seed->the_ptr()->layer() ) {
T getParameter(std::string const &) const
const math::XYZPoint & position() const
cluster centroid position
CaloCluster_iterator preshowerClustersBegin() const
fist iterator over PreshowerCluster constituents
std::auto_ptr< reco::SuperClusterCollection > superClustersEB_
double etawidthSuperClusterBarrel_
PFECALSuperClusterAlgo()
constructor
double phiwidthSuperClusterEndcap_
double correctedEnergy() const
void addHitAndFraction(DetId id, float fraction)
edm::EDGetTokenT< edm::View< reco::PFCluster > > inputTagPFClusters_
bool inDynamicDPhiWindow(const float seedEta, const float seedPhi, const float ClustE, const float ClusEta, const float clusPhi)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
void setPreshowerEnergyPlane2(double preshowerEnergy2)
T const * get() const
Returns C++ pointer to the item.
double threshPFClusterSeedBarrel_
bool applyCrackCorrections_
double pflowPhiWidth() const
const reco::BeamSpot * beamSpot_
void setEnergy(double energy)
const reco::PFCluster::EEtoPSAssociation * EEtoPS_
double Phi_mpi_pi(double x)
void setSeed(const CaloClusterPtr &r)
list of used xtals by DetId // now inherited by CaloCluster
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
void setPhiWidth(double pw)
double etawidthSuperClusterEndcap_
double threshPFClusterSeedEndcap_
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)
std::vector< CalibratedClusterPtr > CalibratedClusterPtrVector
std::shared_ptr< CalibratedPFCluster > CalibratedClusterPtr
MVATrainerComputer * calib
void buildAllSuperClusters(CalibratedClusterPtrVector &, double seedthresh)
std::shared_ptr< PFEnergyCalibration > _pfEnergyCalibration
void setTokens(const edm::ParameterSet &, edm::ConsumesCollector &&)
double dPhi(double phi1, double phi2)
std::vector< SuperCluster > SuperClusterCollection
collection of SuperCluser objectr
bool doSatelliteClusterMerge_
void setCorrectedEnergy(double cenergy)
std::unique_ptr< PFSCRegressionCalc > regr_
std::vector< std::pair< CaloClusterPtr::key_type, edm::Ptr< PFCluster > > > EEtoPSAssociation
Abs< T >::type abs(const T &t)
CalibratedClusterPtrVector _clustersEE
double threshSuperClusterEt_
double energy() const
cluster energy
double satelliteThreshold_
double fractionForMajority_
double threshPFClusterBarrel_
double phiwidthSuperClusterBarrel_
double rawEnergy() const
raw uncorrected energy (sum of energies of component BasicClusters)
edm::EDGetTokenT< reco::PFCluster::EEtoPSAssociation > inputTagPFClustersES_
T const * product() const
edm::EDGetTokenT< reco::BeamSpot > inputTagBeamSpot_
XYZPointD XYZPoint
point in space with cartesian internal representation
void addPreshowerCluster(const CaloClusterPtr &r)
add reference to constituent BasicCluster
std::auto_ptr< reco::SuperClusterCollection > superClustersEE_
clustering_type _clustype
double threshPFClusterEndcap_
void buildSuperCluster(CalibratedClusterPtr &, CalibratedClusterPtrVector &)
bool GreaterByE(const T &a1, const T &a2)
void addCluster(const CaloClusterPtr &r)
add reference to constituent BasicCluster
static int position[264][3]
const Point & position() const
position
void setPreshowerEnergyPlane1(double preshowerEnergy1)
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)
CaloCluster_iterator preshowerClustersEnd() const
last iterator over PreshowerCluster constituents
SCRegressionCalculator< BaselinePFSCRegression > PFSCRegressionCalc
void setPreshowerEnergy(double preshowerEnergy)