CMS 3D CMS Logo

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

#include <PFClusterShapeAlgo.h>

Public Member Functions

reco::ClusterShapeCollectionmakeClusterShapes (edm::Handle< reco::PFClusterCollection > clusterHandle, edm::Handle< reco::PFRecHitCollection > rechitHandle, const CaloSubdetectorGeometry *barrelGeo_p, const CaloSubdetectorTopology *barrelTop_p, const CaloSubdetectorGeometry *endcapGeo_p, const CaloSubdetectorTopology *endcapTop_p)
 
 PFClusterShapeAlgo (bool useFractions, double w0)
 
 ~PFClusterShapeAlgo ()
 

Private Types

enum  Direction {
  N, NE, E, SE,
  S, SW, W, NW
}
 
enum  GeomTopoIndex { BARREL = 0, ENDCAP = 1 }
 
typedef std::map< DetId, RecHitWithFractionRecHitMap
 

Private Member Functions

double addMapEnergies (int etaIndexLow, int etaIndexHigh, int phiIndexLow, int phiIndexHigh)
 
void covariances ()
 
void fill5x5Map ()
 
void find_e2x2 ()
 
void find_e2x5Bottom ()
 
void find_e2x5Left ()
 
void find_e2x5Right ()
 
void find_e2x5Top ()
 
void find_e3x2 ()
 
void find_e3x3 ()
 
void find_e4x4 ()
 
void find_e5x5 ()
 
void find_eMax_e2nd ()
 
int findPFRHIndexFromDetId (unsigned int id)
 
const reco::PFRecHitFractiongetFractionFromDetId (const DetId &id)
 
reco::ClusterShape makeClusterShape ()
 

Private Attributes

double covEtaEta_
 
double covEtaPhi_
 
double covPhiPhi_
 
reco::PFClusterRef currentCluster_p
 
unsigned int currentClusterIndex_
 
edm::Handle< reco::PFRecHitCollectioncurrentRecHit_v_p
 
double e2nd_
 
DetId e2ndId_
 
double e2x2_
 
double e2x5Bottom_
 
double e2x5Left_
 
double e2x5Right_
 
double e2x5Top_
 
double e3x2_
 
double e3x2Ratio_
 
double e3x3_
 
double e4x4_
 
double e5x5_
 
double eMax_
 
Direction eMaxDir
 
DetId eMaxId_
 
unsigned int geomIndex
 
std::vector< const CaloSubdetectorGeometry * > geomVector
 
RecHitWithFraction map5x5 [5][5]
 
math::XYZVector meanPosition_
 
unsigned int topoIndex
 
std::vector< const CaloSubdetectorTopology * > topoVector
 
double totalE_
 
bool useFractions_
 
double w0_
 

Detailed Description

Definition at line 37 of file PFClusterShapeAlgo.h.

Member Typedef Documentation

Definition at line 39 of file PFClusterShapeAlgo.h.

Member Enumeration Documentation

Enumerator
BARREL 
ENDCAP 

Definition at line 42 of file PFClusterShapeAlgo.h.

Constructor & Destructor Documentation

PFClusterShapeAlgo::PFClusterShapeAlgo ( bool  useFractions,
double  w0 
)
explicit
PFClusterShapeAlgo::~PFClusterShapeAlgo ( )

Definition at line 9 of file PFClusterShapeAlgo.cc.

10 {
11 }

Member Function Documentation

double PFClusterShapeAlgo::addMapEnergies ( int  etaIndexLow,
int  etaIndexHigh,
int  phiIndexLow,
int  phiIndexHigh 
)
private

Definition at line 208 of file PFClusterShapeAlgo.cc.

References RecHitWithFraction::energy, mps_fire::i, and map5x5.

Referenced by find_e2x2(), find_e2x5Bottom(), find_e2x5Left(), find_e2x5Right(), find_e2x5Top(), find_e3x2(), find_e3x3(), find_e4x4(), and find_e5x5().

209 {
210  const int etaLow = etaIndexLow + 2;
211  const int etaHigh = etaIndexHigh + 2;
212  const int phiLow = phiIndexLow + 2;
213  const int phiHigh = phiIndexHigh + 2;
214 
215  double energy = 0;
216 
217  for (int i = etaLow; i <= etaHigh; ++i)
218  {
219  for (int j = phiLow; j <= phiHigh; ++j)
220  {
221  energy += map5x5[i][j].energy;
222  }
223  }
224  return energy;
225 }
RecHitWithFraction map5x5[5][5]
void PFClusterShapeAlgo::covariances ( )
private

Definition at line 299 of file PFClusterShapeAlgo.cc.

References covEtaEta_, covEtaPhi_, covPhiPhi_, pfDeepCMVADiscriminatorsJetTags_cfi::denominator, particleFlow_cfi::dEta, particleFlow_cfi::dPhi, mps_fire::i, cmsBatch::log, map5x5, SiStripPI::max, meanPosition_, Geom::pi(), position, totalE_, Geom::twoPi(), w, and w0_.

Referenced by makeClusterShape().

300 {
301  double numeratorEtaEta = 0;
302  double numeratorEtaPhi = 0;
303  double numeratorPhiPhi = 0;
304  double denominator = 0;
305 
306  for (int i = 0; i < 5; ++i)
307  {
308  for (int j = 0; j < 5; ++j)
309  {
310  const math::XYZVector & crystalPosition(map5x5[i][j].position);
311 
312  double dPhi = crystalPosition.phi() - meanPosition_.phi();
313  if (dPhi > + Geom::pi()) { dPhi = Geom::twoPi() - dPhi; }
314  if (dPhi < - Geom::pi()) { dPhi = Geom::twoPi() - dPhi; }
315 
316  const double dEta = crystalPosition.eta() - meanPosition_.eta();
317 
318  const double w = std::max(0.0, w0_ + log(map5x5[i][j].energy / totalE_));
319 
320  denominator += w;
321  numeratorEtaEta += w * dEta * dEta;
322  numeratorEtaPhi += w * dEta * dPhi;
323  numeratorPhiPhi += w * dPhi * dPhi;
324  }
325  }
326 
327  covEtaEta_ = numeratorEtaEta / denominator;
328  covEtaPhi_ = numeratorEtaPhi / denominator;
329  covPhiPhi_ = numeratorPhiPhi / denominator;
330 }
const double w
Definition: UKUtility.cc:23
math::XYZVector meanPosition_
RecHitWithFraction map5x5[5][5]
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
static int position[264][3]
Definition: ReadPGInfo.cc:509
constexpr double pi()
Definition: Pi.h:31
constexpr double twoPi()
Definition: Pi.h:32
void PFClusterShapeAlgo::fill5x5Map ( )
private

Definition at line 155 of file PFClusterShapeAlgo.cc.

References currentRecHit_v_p, RecHitWithFraction::detId, eMaxId_, RecHitWithFraction::energy, findPFRHIndexFromDetId(), dedxEstimators_cff::fraction, reco::PFRecHitFraction::fraction(), getFractionFromDetId(), CaloNavigator< T, TOPO >::home(), mps_fire::i, map5x5, meanPosition_, CaloNavigator< T, TOPO >::offsetBy(), RecHitWithFraction::position, position, topoIndex, topoVector, totalE_, and useFractions_.

Referenced by makeClusterShape().

156 {
157  // first get a navigator to the central element
159 
160  meanPosition_ = math::XYZVector(0.0, 0.0, 0.0);
161  totalE_ = 0;
162 
163  for (int i = 0; i < 5; ++i)
164  {
165  for (int j = 0; j < 5; ++j)
166  {
167  position.home();
168  position.offsetBy(i - 2, j - 2);
169 
170  RecHitWithFraction newEntry;
171  newEntry.detId = DetId(0);
172  newEntry.energy = 0.0;
173  newEntry.position = math::XYZVector(0, 0, 0);
174 
175  if (*position != DetId(0)) // if this is a valid detId...
176  {
177  // ...find the corresponding PFRecHit index
178  const int index = findPFRHIndexFromDetId((*position).rawId());
179 
180  if (index >= 0) // if a PFRecHit exists for this detId
181  {
182  double fraction = 1.0;
183  if (useFractions_) // if the algorithm should use fractions
184  {
185  fraction = 0.0;
186  const reco::PFRecHitFraction * fraction_p = getFractionFromDetId(*position);
187  if (fraction_p) { fraction = fraction_p->fraction(); }
188  }
189 
190  const reco::PFRecHitRef rhRef(currentRecHit_v_p, index);
191  const math::XYZVector crystalPosition(rhRef->position());
192  const double energyFraction = rhRef->energy() * fraction;
193 
194  newEntry.detId = *position;
195  newEntry.energy = energyFraction;
196  newEntry.position = crystalPosition;
197 
198  meanPosition_ = meanPosition_ + crystalPosition * energyFraction;
199  totalE_ += energyFraction;
200  }
201  }
202  map5x5[i][j] = newEntry;
203  }
204  }
205  meanPosition_ /= totalE_;
206 }
math::XYZVector position
Fraction of a PFRecHit (rechits can be shared between several PFCluster&#39;s)
double fraction() const
int findPFRHIndexFromDetId(unsigned int id)
math::XYZVector meanPosition_
T offsetBy(int deltaX, int deltaY) const
Free movement of arbitray steps.
Definition: CaloNavigator.h:80
RecHitWithFraction map5x5[5][5]
void home() const
move the navigator back to the starting point
Definition: DetId.h:18
edm::Handle< reco::PFRecHitCollection > currentRecHit_v_p
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
std::vector< const CaloSubdetectorTopology * > topoVector
static int position[264][3]
Definition: ReadPGInfo.cc:509
const reco::PFRecHitFraction * getFractionFromDetId(const DetId &id)
void PFClusterShapeAlgo::find_e2x2 ( )
private

Definition at line 242 of file PFClusterShapeAlgo.cc.

References addMapEnergies(), e2x2_, eMaxDir, NE, NW, SE, and SW.

Referenced by makeClusterShape().

243 {
244  std::map<double, Direction> directionMap;
245 
246  directionMap[addMapEnergies(-1, +0, -1, +0)] = SE;
247  directionMap[addMapEnergies(-1, +0, +0, +1)] = NE;
248  directionMap[addMapEnergies(+0, +1, -1, +0)] = SW;
249  directionMap[addMapEnergies(+0, +1, +0, +1)] = NW;
250 
251  const std::map<double, Direction>::reverse_iterator eMaxDir_it = directionMap.rbegin();
252 
253  eMaxDir = eMaxDir_it->second;
254 
255  e2x2_ = eMaxDir_it->first;
256 }
double addMapEnergies(int etaIndexLow, int etaIndexHigh, int phiIndexLow, int phiIndexHigh)
void PFClusterShapeAlgo::find_e2x5Bottom ( )
private

Definition at line 232 of file PFClusterShapeAlgo.cc.

References addMapEnergies(), and e2x5Bottom_.

Referenced by makeClusterShape().

232 { e2x5Bottom_ = addMapEnergies(+1, +2, -2, +2); }
double addMapEnergies(int etaIndexLow, int etaIndexHigh, int phiIndexLow, int phiIndexHigh)
void PFClusterShapeAlgo::find_e2x5Left ( )
private

Definition at line 230 of file PFClusterShapeAlgo.cc.

References addMapEnergies(), and e2x5Left_.

Referenced by makeClusterShape().

230 { e2x5Left_ = addMapEnergies(-2, +2, -2, -1); }
double addMapEnergies(int etaIndexLow, int etaIndexHigh, int phiIndexLow, int phiIndexHigh)
void PFClusterShapeAlgo::find_e2x5Right ( )
private

Definition at line 229 of file PFClusterShapeAlgo.cc.

References addMapEnergies(), and e2x5Right_.

Referenced by makeClusterShape().

229 { e2x5Right_ = addMapEnergies(-2, +2, +1, +2); }
double addMapEnergies(int etaIndexLow, int etaIndexHigh, int phiIndexLow, int phiIndexHigh)
void PFClusterShapeAlgo::find_e2x5Top ( )
private

Definition at line 231 of file PFClusterShapeAlgo.cc.

References addMapEnergies(), and e2x5Top_.

Referenced by makeClusterShape().

231 { e2x5Top_ = addMapEnergies(-2, -1, -2, +2); }
double addMapEnergies(int etaIndexLow, int etaIndexHigh, int phiIndexLow, int phiIndexHigh)
void PFClusterShapeAlgo::find_e3x2 ( )
private

Definition at line 258 of file PFClusterShapeAlgo.cc.

References addMapEnergies(), pfDeepCMVADiscriminatorsJetTags_cfi::denominator, dir, E, e3x2_, e3x2Ratio_, RecHitWithFraction::energy, map5x5, N, pfDeepCMVADiscriminatorsJetTags_cfi::numerator, S, and W.

Referenced by makeClusterShape().

259 {
260  // Find the direction of the highest energy neighbour
261  std::map<double, Direction> directionMap;
262  directionMap[map5x5[2][3].energy] = N;
263  directionMap[map5x5[2][1].energy] = S;
264  directionMap[map5x5[1][2].energy] = E;
265  directionMap[map5x5[3][2].energy] = W;
266  // Maps are sorted in ascending order - get the last element
267  const Direction dir = directionMap.rbegin()->second;
268 
269  if (dir == N)
270  {
271  e3x2_ = addMapEnergies(-1, +1, +0, +1);
272  const double numerator = map5x5[3][2].energy + map5x5[1][2].energy;
273  const double denominator = map5x5[1][3].energy + map5x5[3][3].energy + 0.5;
274  e3x2Ratio_ = numerator / denominator;
275  }
276  else if (dir == S)
277  {
278  e3x2_ = addMapEnergies(-1, +1, -1, +0);
279  const double numerator = map5x5[3][2].energy + map5x5[1][2].energy;
280  const double denominator = map5x5[1][1].energy + map5x5[3][1].energy + 0.5;
281  e3x2Ratio_ = numerator / denominator;
282  }
283  else if (dir == W)
284  {
285  e3x2_ = addMapEnergies(+0, +1, -1, +1);
286  const double numerator = map5x5[2][3].energy + map5x5[2][1].energy;
287  const double denominator = map5x5[3][3].energy + map5x5[3][1].energy + 0.5;
288  e3x2Ratio_ = numerator / denominator;
289  }
290  else if (dir == E)
291  {
292  e3x2_ = addMapEnergies(-1, +0, -1, +1);
293  const double numerator = map5x5[2][3].energy + map5x5[2][1].energy;
294  const double denominator = map5x5[1][1].energy + map5x5[1][3].energy + 0.5;
295  e3x2Ratio_ = numerator / denominator;
296  }
297 }
RecHitWithFraction map5x5[5][5]
double addMapEnergies(int etaIndexLow, int etaIndexHigh, int phiIndexLow, int phiIndexHigh)
dbl *** dir
Definition: mlp_gen.cc:35
void PFClusterShapeAlgo::find_e3x3 ( )
private

Definition at line 227 of file PFClusterShapeAlgo.cc.

References addMapEnergies(), and e3x3_.

Referenced by makeClusterShape().

227 { e3x3_ = addMapEnergies(-1, +1, -1, +1); }
double addMapEnergies(int etaIndexLow, int etaIndexHigh, int phiIndexLow, int phiIndexHigh)
void PFClusterShapeAlgo::find_e4x4 ( )
private

Definition at line 234 of file PFClusterShapeAlgo.cc.

References addMapEnergies(), e4x4_, eMaxDir, NE, NW, SE, and SW.

Referenced by makeClusterShape().

235 {
236  if (eMaxDir == SE) { e4x4_ = addMapEnergies(-2, +1, -2, +1); return; }
237  if (eMaxDir == NE) { e4x4_ = addMapEnergies(-2, +1, -1, +2); return; }
238  if (eMaxDir == SW) { e4x4_ = addMapEnergies(-1, +2, -2, +1); return; }
239  if (eMaxDir == NW) { e4x4_ = addMapEnergies(-1, +2, -1, +2); return; }
240 }
double addMapEnergies(int etaIndexLow, int etaIndexHigh, int phiIndexLow, int phiIndexHigh)
void PFClusterShapeAlgo::find_e5x5 ( )
private

Definition at line 228 of file PFClusterShapeAlgo.cc.

References addMapEnergies(), and e5x5_.

Referenced by makeClusterShape().

228 { e5x5_ = addMapEnergies(-2, +2, -2, +2); }
double addMapEnergies(int etaIndexLow, int etaIndexHigh, int phiIndexLow, int phiIndexHigh)
void PFClusterShapeAlgo::find_eMax_e2nd ( )
private

Definition at line 96 of file PFClusterShapeAlgo.cc.

References currentCluster_p, e2nd_, e2ndId_, eMax_, and eMaxId_.

Referenced by makeClusterShape().

97 {
98  std::map<double, DetId> energyMap;
99 
100  // First get the RecHitFractions:
101  const std::vector<reco::PFRecHitFraction> & fraction_v = currentCluster_p->recHitFractions();
102  // For every one of them...
103  for (std::vector<reco::PFRecHitFraction>::const_iterator it = fraction_v.begin(); it != fraction_v.end(); ++it)
104  {
105  // ...find the corresponding rechit:
106  // const reco::PFRecHit & rechit = (*currentRecHit_v_p)[it->recHitIndex()];
107  const reco::PFRecHitRef rechit = it->recHitRef();
108  // ...and DetId:
109  const DetId rechitDetId = DetId(rechit->detId());
110  // Make the new Pair and put it in the map:
111  energyMap[rechit->energy()] = rechitDetId;
112  }
113  // maps are sorted in ascending order so get the last two elements:
114  std::map<double, DetId>::reverse_iterator it = energyMap.rbegin();
115  eMax_ = it->first;
116  eMaxId_ = it->second;
117  it++;
118  e2nd_ = it->first;
119  e2ndId_ = it->second;
120 }
reco::PFClusterRef currentCluster_p
Definition: DetId.h:18
int PFClusterShapeAlgo::findPFRHIndexFromDetId ( unsigned int  id)
private

Definition at line 122 of file PFClusterShapeAlgo.cc.

References currentRecHit_v_p, triggerObjects_cff::id, and gen::k.

Referenced by fill5x5Map().

123 {
124  int index = -1; // need some negative number
125  for (unsigned int k = 0; k < currentRecHit_v_p->size(); ++k)
126  {
127  if ((*currentRecHit_v_p)[k].detId() == id)
128  {
129  index = static_cast<int>(k);
130  break;
131  }
132  }
133  return index;
134 }
int k[5][pyjets_maxn]
edm::Handle< reco::PFRecHitCollection > currentRecHit_v_p
const reco::PFRecHitFraction * PFClusterShapeAlgo::getFractionFromDetId ( const DetId id)
private

Definition at line 137 of file PFClusterShapeAlgo.cc.

References currentCluster_p.

Referenced by fill5x5Map().

138 {
139  const std::vector< reco::PFRecHitFraction > & fraction_v = currentCluster_p->recHitFractions();
140  for (std::vector<reco::PFRecHitFraction>::const_iterator it = fraction_v.begin(); it != fraction_v.end(); ++it)
141  {
142  //const unsigned int rhIndex = it->recHitIndex();
143  //reco::PFRecHitRef rh_p(currentRecHit_v_p, rhIndex);
144  const reco::PFRecHitRef rh_p = it->recHitRef();
145  const DetId rhDetId = DetId(rh_p->detId());
146  if (rhDetId == id)
147  {
148  return &(*it);
149  }
150  }
151  return nullptr;
152 }
reco::PFClusterRef currentCluster_p
Definition: DetId.h:18
reco::ClusterShape PFClusterShapeAlgo::makeClusterShape ( )
private

Definition at line 60 of file PFClusterShapeAlgo.cc.

References covariances(), covEtaEta_, covEtaPhi_, covPhiPhi_, e2nd_, e2ndId_, e2x2_, e2x5Bottom_, e2x5Left_, e2x5Right_, e2x5Top_, e3x2_, e3x2Ratio_, e3x3_, e4x4_, e5x5_, eMax_, eMaxId_, fill5x5Map(), find_e2x2(), find_e2x5Bottom(), find_e2x5Left(), find_e2x5Right(), find_e2x5Top(), find_e3x2(), find_e3x3(), find_e4x4(), find_e5x5(), and find_eMax_e2nd().

Referenced by makeClusterShapes().

61 {
63  fill5x5Map();
64 
65  find_e2x2();
66  find_e3x2();
67  find_e3x3();
68  find_e4x4();
69  find_e5x5();
70 
72  find_e2x5Left();
73  find_e2x5Top();
75 
76  covariances();
77 
78  double dummyLAT = 0;
79  double dummyEtaLAT = 0;
80  double dummyPhiLAT = 0;
81  double dummyA20 = 0;
82  double dummyA42 = 0;
83 
84  std::vector<double> dummyEnergyBasketFractionEta_v;
85  std::vector<double> dummyEnergyBasketFractionPhi_v;
86 
91  dummyLAT, dummyEtaLAT, dummyPhiLAT, dummyA20, dummyA42,
92  dummyEnergyBasketFractionEta_v, dummyEnergyBasketFractionPhi_v);
93 }
reco::ClusterShapeCollection * PFClusterShapeAlgo::makeClusterShapes ( edm::Handle< reco::PFClusterCollection clusterHandle,
edm::Handle< reco::PFRecHitCollection rechitHandle,
const CaloSubdetectorGeometry barrelGeo_p,
const CaloSubdetectorTopology barrelTop_p,
const CaloSubdetectorGeometry endcapGeo_p,
const CaloSubdetectorTopology endcapTop_p 
)

Definition at line 14 of file PFClusterShapeAlgo.cc.

References BARREL, currentCluster_p, currentClusterIndex_, currentRecHit_v_p, ENDCAP, geomIndex, geomVector, mps_fire::i, makeClusterShape(), topoIndex, and topoVector.

20 {
21  static const float etaEndOfBarrel = 1.497;
22 
23  topoVector.push_back(the_barrelTop_p);
24  topoVector.push_back(the_endcapTop_p);
25  geomVector.push_back(the_barrelGeo_p);
26  geomVector.push_back(the_endcapGeo_p);
27 
29 
30  currentRecHit_v_p = rechitHandle;
31 
32  for (unsigned int i = 0; i < clusterHandle->size(); ++i)
33  {
34  // Make each cluster the "current" cluster
35  currentCluster_p = reco::PFClusterRef(clusterHandle, i);
37 
38  // Find the right topology to use with this cluster
39  topoIndex = BARREL;
40  geomIndex = BARREL;
41  const math::XYZVector currentClusterPos(currentCluster_p->position());
42  if (fabs(currentClusterPos.eta()) > etaEndOfBarrel)
43  {
44  topoIndex = ENDCAP;
45  geomIndex = ENDCAP;
46  }
47 
48  // Create the clustershape and push it into the vector
49  shape_v_p->push_back(makeClusterShape());
50  }
51 
52  topoVector.clear();
53  topoVector.clear();
54  geomVector.clear();
55  geomVector.clear();
56 
57  return shape_v_p;
58 }
reco::ClusterShape makeClusterShape()
std::vector< ClusterShape > ClusterShapeCollection
collection of ClusterShape objects
std::vector< const CaloSubdetectorGeometry * > geomVector
edm::Ref< PFClusterCollection > PFClusterRef
persistent reference to PFCluster objects
Definition: PFClusterFwd.h:15
reco::PFClusterRef currentCluster_p
unsigned int currentClusterIndex_
edm::Handle< reco::PFRecHitCollection > currentRecHit_v_p
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
std::vector< const CaloSubdetectorTopology * > topoVector

Member Data Documentation

double PFClusterShapeAlgo::covEtaEta_
private

Definition at line 81 of file PFClusterShapeAlgo.h.

Referenced by covariances(), and makeClusterShape().

double PFClusterShapeAlgo::covEtaPhi_
private

Definition at line 81 of file PFClusterShapeAlgo.h.

Referenced by covariances(), and makeClusterShape().

double PFClusterShapeAlgo::covPhiPhi_
private

Definition at line 81 of file PFClusterShapeAlgo.h.

Referenced by covariances(), and makeClusterShape().

reco::PFClusterRef PFClusterShapeAlgo::currentCluster_p
private

Definition at line 64 of file PFClusterShapeAlgo.h.

Referenced by find_eMax_e2nd(), getFractionFromDetId(), and makeClusterShapes().

unsigned int PFClusterShapeAlgo::currentClusterIndex_
private

Definition at line 63 of file PFClusterShapeAlgo.h.

Referenced by makeClusterShapes().

edm::Handle<reco::PFRecHitCollection> PFClusterShapeAlgo::currentRecHit_v_p
private

Definition at line 65 of file PFClusterShapeAlgo.h.

Referenced by fill5x5Map(), findPFRHIndexFromDetId(), and makeClusterShapes().

double PFClusterShapeAlgo::e2nd_
private

Definition at line 79 of file PFClusterShapeAlgo.h.

Referenced by find_eMax_e2nd(), and makeClusterShape().

DetId PFClusterShapeAlgo::e2ndId_
private

Definition at line 78 of file PFClusterShapeAlgo.h.

Referenced by find_eMax_e2nd(), and makeClusterShape().

double PFClusterShapeAlgo::e2x2_
private

Definition at line 80 of file PFClusterShapeAlgo.h.

Referenced by find_e2x2(), and makeClusterShape().

double PFClusterShapeAlgo::e2x5Bottom_
private

Definition at line 80 of file PFClusterShapeAlgo.h.

Referenced by find_e2x5Bottom(), and makeClusterShape().

double PFClusterShapeAlgo::e2x5Left_
private

Definition at line 80 of file PFClusterShapeAlgo.h.

Referenced by find_e2x5Left(), and makeClusterShape().

double PFClusterShapeAlgo::e2x5Right_
private

Definition at line 80 of file PFClusterShapeAlgo.h.

Referenced by find_e2x5Right(), and makeClusterShape().

double PFClusterShapeAlgo::e2x5Top_
private

Definition at line 80 of file PFClusterShapeAlgo.h.

Referenced by find_e2x5Top(), and makeClusterShape().

double PFClusterShapeAlgo::e3x2_
private

Definition at line 80 of file PFClusterShapeAlgo.h.

Referenced by find_e3x2(), and makeClusterShape().

double PFClusterShapeAlgo::e3x2Ratio_
private

Definition at line 80 of file PFClusterShapeAlgo.h.

Referenced by find_e3x2(), and makeClusterShape().

double PFClusterShapeAlgo::e3x3_
private

Definition at line 80 of file PFClusterShapeAlgo.h.

Referenced by find_e3x3(), and makeClusterShape().

double PFClusterShapeAlgo::e4x4_
private

Definition at line 80 of file PFClusterShapeAlgo.h.

Referenced by find_e4x4(), and makeClusterShape().

double PFClusterShapeAlgo::e5x5_
private

Definition at line 80 of file PFClusterShapeAlgo.h.

Referenced by find_e5x5(), and makeClusterShape().

double PFClusterShapeAlgo::eMax_
private

Definition at line 79 of file PFClusterShapeAlgo.h.

Referenced by find_eMax_e2nd(), and makeClusterShape().

Direction PFClusterShapeAlgo::eMaxDir
private

Definition at line 76 of file PFClusterShapeAlgo.h.

Referenced by find_e2x2(), and find_e4x4().

DetId PFClusterShapeAlgo::eMaxId_
private

Definition at line 78 of file PFClusterShapeAlgo.h.

Referenced by fill5x5Map(), find_eMax_e2nd(), and makeClusterShape().

unsigned int PFClusterShapeAlgo::geomIndex
private

Definition at line 69 of file PFClusterShapeAlgo.h.

Referenced by makeClusterShapes().

std::vector<const CaloSubdetectorGeometry *> PFClusterShapeAlgo::geomVector
private

Definition at line 70 of file PFClusterShapeAlgo.h.

Referenced by makeClusterShapes().

RecHitWithFraction PFClusterShapeAlgo::map5x5[5][5]
private

Definition at line 72 of file PFClusterShapeAlgo.h.

Referenced by addMapEnergies(), covariances(), fill5x5Map(), and find_e3x2().

math::XYZVector PFClusterShapeAlgo::meanPosition_
private

Definition at line 73 of file PFClusterShapeAlgo.h.

Referenced by covariances(), and fill5x5Map().

unsigned int PFClusterShapeAlgo::topoIndex
private

Definition at line 67 of file PFClusterShapeAlgo.h.

Referenced by fill5x5Map(), and makeClusterShapes().

std::vector<const CaloSubdetectorTopology *> PFClusterShapeAlgo::topoVector
private

Definition at line 68 of file PFClusterShapeAlgo.h.

Referenced by fill5x5Map(), and makeClusterShapes().

double PFClusterShapeAlgo::totalE_
private

Definition at line 74 of file PFClusterShapeAlgo.h.

Referenced by covariances(), and fill5x5Map().

bool PFClusterShapeAlgo::useFractions_
private

Definition at line 60 of file PFClusterShapeAlgo.h.

Referenced by fill5x5Map(), and PFClusterShapeAlgo().

double PFClusterShapeAlgo::w0_
private

Definition at line 61 of file PFClusterShapeAlgo.h.

Referenced by covariances(), and PFClusterShapeAlgo().