CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
Multi5x5BremRecoveryClusterAlgo Class Reference

#include <Multi5x5BremRecoveryClusterAlgo.h>

Public Member Functions

reco::SuperClusterCollection makeSuperClusters (const reco::CaloClusterPtrVector &clusters)
 
 Multi5x5BremRecoveryClusterAlgo (const edm::ParameterSet &bremRecoveryPset, double eb_sc_road_etasize=0.06, double eb_sc_road_phisize=0.80, double ec_sc_road_etasize=0.14, double ec_sc_road_phisize=0.40, bool dynamicPhiRoad=true, double theSeedTransverseEnergyThreshold=0.40)
 
 ~Multi5x5BremRecoveryClusterAlgo ()
 

Private Member Functions

void makeIslandSuperClusters (reco::CaloClusterPtrVector &clusters_v, double etaRoad, double phiRoad)
 

Private Attributes

bool dynamicPhiRoad_
 
double eb_rdeta_
 
double eb_rdphi_
 
double ec_rdeta_
 
double ec_rdphi_
 
BremRecoveryPhiRoadAlgophiRoadAlgo_
 
double seedTransverseEnergyThreshold
 
reco::SuperClusterCollection superclusters_v
 

Detailed Description

Definition at line 24 of file Multi5x5BremRecoveryClusterAlgo.h.

Constructor & Destructor Documentation

Multi5x5BremRecoveryClusterAlgo::Multi5x5BremRecoveryClusterAlgo ( const edm::ParameterSet bremRecoveryPset,
double  eb_sc_road_etasize = 0.06,
double  eb_sc_road_phisize = 0.80,
double  ec_sc_road_etasize = 0.14,
double  ec_sc_road_phisize = 0.40,
bool  dynamicPhiRoad = true,
double  theSeedTransverseEnergyThreshold = 0.40 
)
inline

Definition at line 29 of file Multi5x5BremRecoveryClusterAlgo.h.

References cosmicSuperClusters_cfi::dynamicPhiRoad, dynamicPhiRoad_, eb_rdeta_, eb_rdphi_, ec_rdeta_, ec_rdphi_, phiRoadAlgo_, and seedTransverseEnergyThreshold.

37  {
38  // e*_rdeta_ and e*_rdphi_ are half the total window
39  // because they correspond to one direction (positive or negative)
40  eb_rdeta_ = eb_sc_road_etasize / 2;
41  eb_rdphi_ = eb_sc_road_phisize / 2;
42  ec_rdeta_ = ec_sc_road_etasize / 2;
43  ec_rdphi_ = ec_sc_road_phisize / 2;
44 
45  seedTransverseEnergyThreshold = theSeedTransverseEnergyThreshold;
47  if (dynamicPhiRoad_) phiRoadAlgo_ = new BremRecoveryPhiRoadAlgo(bremRecoveryPset);
48 
49  }
Multi5x5BremRecoveryClusterAlgo::~Multi5x5BremRecoveryClusterAlgo ( )
inline

Member Function Documentation

void Multi5x5BremRecoveryClusterAlgo::makeIslandSuperClusters ( reco::CaloClusterPtrVector clusters_v,
double  etaRoad,
double  phiRoad 
)
private

Definition at line 41 of file Multi5x5BremRecoveryClusterAlgo.cc.

References funct::abs(), edm::PtrVectorBase::clear(), dynamicPhiRoad_, edm::PtrVectorBase::empty(), BremRecoveryPhiRoadAlgo::endcapPhiRoad(), reco::CaloCluster::energy(), stringResolutionProvider_cfi::et, PVValHelper::eta, objects.autophobj::float, mps_fire::i, LogTrace, M_PI, match(), p1, p2, phi, phiRoadAlgo_, reco::CaloCluster::position(), position, edm::PtrVector< T >::push_back(), seedTransverseEnergyThreshold, funct::sin(), edm::PtrVectorBase::size(), superclusters_v, theta(), mitigatedMETSequence_cff::U, X, DOFs::Y, and DOFs::Z.

Referenced by makeSuperClusters(), and ~Multi5x5BremRecoveryClusterAlgo().

43 {
44  if(clusters_v.empty()) return;
45 
46  const auto clustersSize = clusters_v.size();
47  assert(clustersSize>0);
48 
49  bool usedSeed[clustersSize];
50  for (auto ic=0U; ic<clustersSize; ++ic) usedSeed[ic]=false;
51 
52  float eta[clustersSize], phi[clustersSize], et[clustersSize];
53  for (auto ic=0U; ic<clustersSize; ++ic) {
54  eta[ic]=clusters_v[ic]->eta();
55  phi[ic]=clusters_v[ic]->phi();
56  et[ic]=clusters_v[ic]->energy() * sin(clusters_v[ic]->position().theta());
57  }
58 
59  for (auto is=0U; is<clustersSize; ++is) {
60  // check this seed was not already used
61  if (usedSeed[is]) continue;
62  auto const & currentSeed = clusters_v[is];
63 
64  // Does our highest energy cluster have high enough energy?
65  // changed this to continue from break (to be robust against the order of sorting of the seed clusters)
66  if (et[is] < seedTransverseEnergyThreshold) continue;
67 
68  // if yes, make it a seed for a new SuperCluster:
69  double energy = (currentSeed)->energy();
70  math::XYZVector position_((currentSeed)->position().X(),
71  (currentSeed)->position().Y(),
72  (currentSeed)->position().Z());
73  position_ *= energy;
74  usedSeed[is]=true;
75 
76  LogTrace("EcalClusters") << "*****************************";
77  LogTrace("EcalClusters") << "******NEW SUPERCLUSTER*******";
78  LogTrace("EcalClusters") << "Seed R = " << (currentSeed)->position().Rho();
79 
80  reco::CaloClusterPtrVector constituentClusters;
81  constituentClusters.push_back(currentSeed);
82  auto ic = is + 1;
83 
84  while (ic < clustersSize) {
85  auto const & currentCluster = clusters_v[ic];
86 
87  // if dynamic phi road is enabled then compute the phi road for a cluster
88  // of energy of existing clusters + the candidate cluster.
89  if (dynamicPhiRoad_)
90  phiRoad = phiRoadAlgo_->endcapPhiRoad(energy + (currentCluster)->energy());
91 
92  auto dphi = [](float p1,float p2) {
93  auto dp=std::abs(p1-p2); if (dp>float(M_PI)) dp-=float(2*M_PI);
94  return std::abs(dp);
95  };
96 
97  auto match = [&](int i, int j) {
98  return (dphi(phi[i],phi[j])< phiRoad) & (std::abs(eta[i]-eta[j])< etaRoad);
99  };
100 
101  // does the cluster match the phi road for this candidate supercluster
102  if (!usedSeed[ic] && match(is,ic)) {
103 
104  // add basic cluster to supercluster constituents
105  constituentClusters.push_back(currentCluster);
106  energy += (currentCluster)->energy();
107  position_ += (currentCluster)->energy() * math::XYZVector((currentCluster)->position().X(),
108  (currentCluster)->position().Y(),
109  (currentCluster)->position().Z()
110  );
111  // remove cluster from vector of available clusters
112  usedSeed[ic]=true;
113  LogTrace("EcalClusters") << "Cluster R = " << (currentCluster)->position().Rho();
114  }
115  ++ic;
116  }
117 
118  position_ /= energy;
119 
120  LogTrace("EcalClusters") <<"Final SuperCluster R = " << position_.Rho();
121 
122  reco::SuperCluster newSuperCluster(energy,
123  math::XYZPoint(position_.X(), position_.Y(), position_.Z()),
124  currentSeed,
125  constituentClusters);
126 
127  superclusters_v.push_back(newSuperCluster);
128  LogTrace("EcalClusters") << "created a new supercluster of: ";
129  LogTrace("EcalClusters") << "Energy = " << newSuperCluster.energy();
130  LogTrace("EcalClusters") << "Position in (R, phi, theta) = ("
131  << newSuperCluster.position().Rho() << ", "
132  << newSuperCluster.position().phi() << ", "
133  << newSuperCluster.position().theta() << ")" ;
134 
135 
136  }
137 
138  clusters_v.clear();
139 
140 }
size_type size() const
Size of the RefVector.
Definition: PtrVectorBase.h:74
void push_back(Ptr< T > const &iPtr)
Definition: PtrVector.h:140
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
bool empty() const
Is the RefVector empty.
Definition: PtrVectorBase.h:71
Geom::Theta< T > theta() const
#define X(str)
Definition: MuonsGrabber.cc:48
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double p2[4]
Definition: TauolaWrapper.h:90
#define LogTrace(id)
#define M_PI
reco::SuperClusterCollection superclusters_v
double endcapPhiRoad(double energy)
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
et
define resolution functions of each parameter
void clear()
Clear the PtrVector.
Definition: PtrVectorBase.h:80
double p1[4]
Definition: TauolaWrapper.h:89
static int position[264][3]
Definition: ReadPGInfo.cc:509
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
reco::SuperClusterCollection Multi5x5BremRecoveryClusterAlgo::makeSuperClusters ( const reco::CaloClusterPtrVector clusters)

Definition at line 5 of file Multi5x5BremRecoveryClusterAlgo.cc.

References eb_rdeta_, eb_rdphi_, ec_rdeta_, ec_rdphi_, makeIslandSuperClusters(), reco::CaloCluster::multi5x5, edm::PtrVector< T >::push_back(), and superclusters_v.

Referenced by Multi5x5SuperClusterProducer::produceSuperclustersForECALPart(), and ~Multi5x5BremRecoveryClusterAlgo().

6 {
7 
8  const float etaBorder = 1.479;
9  superclusters_v.clear();
10 
11  // create vectors of references to clusters of a specific origin...
12  reco::CaloClusterPtrVector islandClustersBarrel_v;
13  reco::CaloClusterPtrVector islandClustersEndCap_v;
14 
15  // ...and populate them:
16  for (auto const& cluster_p : clustersCollection)
17  {
18  if (cluster_p->algo() == reco::CaloCluster::multi5x5)
19  {
20  if (fabs(cluster_p->position().eta()) < etaBorder)
21  {
22  islandClustersBarrel_v.push_back(cluster_p);
23  }
24  else
25  {
26  islandClustersEndCap_v.push_back(cluster_p);
27  }
28  }
29  }
30 
31  // make the superclusters from the Barrel clusters - Island
32  makeIslandSuperClusters(islandClustersBarrel_v, eb_rdeta_, eb_rdphi_);
33  // make the superclusters from the EndCap clusters - Island
34  makeIslandSuperClusters(islandClustersEndCap_v, ec_rdeta_, ec_rdphi_);
35 
36  return superclusters_v;
37 }
void makeIslandSuperClusters(reco::CaloClusterPtrVector &clusters_v, double etaRoad, double phiRoad)
void push_back(Ptr< T > const &iPtr)
Definition: PtrVector.h:140
reco::SuperClusterCollection superclusters_v

Member Data Documentation

bool Multi5x5BremRecoveryClusterAlgo::dynamicPhiRoad_
private
double Multi5x5BremRecoveryClusterAlgo::eb_rdeta_
private
double Multi5x5BremRecoveryClusterAlgo::eb_rdphi_
private
double Multi5x5BremRecoveryClusterAlgo::ec_rdeta_
private
double Multi5x5BremRecoveryClusterAlgo::ec_rdphi_
private
BremRecoveryPhiRoadAlgo* Multi5x5BremRecoveryClusterAlgo::phiRoadAlgo_
private
double Multi5x5BremRecoveryClusterAlgo::seedTransverseEnergyThreshold
private
reco::SuperClusterCollection Multi5x5BremRecoveryClusterAlgo::superclusters_v
private

Definition at line 77 of file Multi5x5BremRecoveryClusterAlgo.h.

Referenced by makeIslandSuperClusters(), and makeSuperClusters().