CMS 3D CMS Logo

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

#include <Multi5x5BremRecoveryClusterAlgo.h>

Public Member Functions

reco::SuperClusterCollection makeSuperClusters (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 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;
46  dynamicPhiRoad_ = dynamicPhiRoad;
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 42 of file Multi5x5BremRecoveryClusterAlgo.cc.

References funct::abs(), edm::PtrVectorBase::clear(), reco::dp, dynamicPhiRoad_, BremRecoveryPhiRoadAlgo::endcapPhiRoad(), relval_parameters_module::energy, reco::CaloCluster::energy(), eta(), i, j, 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(), X, and Gflash::Z.

Referenced by makeSuperClusters().

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

Definition at line 5 of file Multi5x5BremRecoveryClusterAlgo.cc.

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

Referenced by Multi5x5SuperClusterProducer::produceSuperclustersForECALPart().

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 (reco::CaloCluster_iterator it = clustersCollection.begin(); it != clustersCollection.end(); it++)
17  {
18  reco::CaloClusterPtr cluster_p = *it;
19  if (cluster_p->algo() == reco::CaloCluster::multi5x5)
20  {
21  if (fabs(cluster_p->position().eta()) < etaBorder)
22  {
23  islandClustersBarrel_v.push_back(cluster_p);
24  }
25  else
26  {
27  islandClustersEndCap_v.push_back(cluster_p);
28  }
29  }
30  }
31 
32  // make the superclusters from the Barrel clusters - Island
33  makeIslandSuperClusters(islandClustersBarrel_v, eb_rdeta_, eb_rdphi_);
34  // make the superclusters from the EndCap clusters - Island
35  makeIslandSuperClusters(islandClustersEndCap_v, ec_rdeta_, ec_rdphi_);
36 
37  return superclusters_v;
38 }
void makeIslandSuperClusters(reco::CaloClusterPtrVector &clusters_v, double etaRoad, double phiRoad)
void push_back(Ptr< T > const &iPtr)
Definition: PtrVector.h:138
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().