CMS 3D CMS Logo

PFClusterShapeAlgo Class Reference

#include <RecoParticleFlow/PFClusterShapeAlgo/interface/PFClusterShapeAlgo.h>

List of all members.

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,
RecHitWithFraction
RecHitMap

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::PFRecHitCollection
currentRecHit_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

typedef std::map<DetId, RecHitWithFraction> PFClusterShapeAlgo::RecHitMap [private]

Definition at line 39 of file PFClusterShapeAlgo.h.


Member Enumeration Documentation

enum PFClusterShapeAlgo::Direction [private]

Enumerator:
N 
NE 
E 
SE 
S 
SW 
W 
NW 

Definition at line 41 of file PFClusterShapeAlgo.h.

00041 { N, NE, E, SE, S, SW, W, NW };

enum PFClusterShapeAlgo::GeomTopoIndex [private]

Enumerator:
BARREL 
ENDCAP 

Definition at line 42 of file PFClusterShapeAlgo.h.

00042 { BARREL = 0, ENDCAP = 1 }; 


Constructor & Destructor Documentation

PFClusterShapeAlgo::PFClusterShapeAlgo ( bool  useFractions,
double  w0 
) [explicit]

Definition at line 3 of file PFClusterShapeAlgo.cc.

References useFractions_, and w0_.

00004 {
00005   useFractions_ = useFractions;
00006   w0_ = w0;
00007 }

PFClusterShapeAlgo::~PFClusterShapeAlgo (  ) 

Definition at line 9 of file PFClusterShapeAlgo.cc.

00010 {
00011 }


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, relval_parameters_module::energy, i, j, and map5x5.

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

00209 {
00210   const int etaLow  = etaIndexLow  + 2;
00211   const int etaHigh = etaIndexHigh + 2;
00212   const int phiLow  = phiIndexLow  + 2;
00213   const int phiHigh = phiIndexHigh + 2;
00214 
00215   double energy = 0;
00216 
00217   for (int i = etaLow; i <= etaHigh; ++i)
00218     {
00219       for (int j = phiLow; j <= phiHigh; ++j)
00220         {
00221           energy += map5x5[i][j].energy;
00222         }
00223     }
00224   return energy;
00225 }

void PFClusterShapeAlgo::covariances (  )  [private]

Definition at line 299 of file PFClusterShapeAlgo.cc.

References covEtaEta_, covEtaPhi_, covPhiPhi_, dPhi(), relval_parameters_module::energy, i, j, funct::log(), map5x5, max, meanPosition_, Geom::pi(), totalE_, Geom::twoPi(), w, and w0_.

Referenced by makeClusterShape().

00300 {
00301   double numeratorEtaEta = 0;
00302   double numeratorEtaPhi = 0;
00303   double numeratorPhiPhi = 0;
00304   double denominator     = 0;
00305 
00306   for (int i = 0; i < 5; ++i)
00307     {
00308       for (int j = 0; j < 5; ++j)
00309         {
00310           const math::XYZVector & crystalPosition(map5x5[i][j].position);
00311           
00312           double dPhi = crystalPosition.phi() - meanPosition_.phi();
00313           if (dPhi > + Geom::pi()) { dPhi = Geom::twoPi() - dPhi; }
00314           if (dPhi < - Geom::pi()) { dPhi = Geom::twoPi() - dPhi; }
00315 
00316           const double dEta = crystalPosition.eta() - meanPosition_.eta();
00317           
00318           const double w = std::max(0.0, w0_ + log(map5x5[i][j].energy / totalE_));
00319           
00320           denominator += w;
00321           numeratorEtaEta += w * dEta * dEta;
00322           numeratorEtaPhi += w * dEta * dPhi;
00323           numeratorPhiPhi += w * dPhi * dPhi;
00324         }
00325     }
00326 
00327   covEtaEta_ = numeratorEtaEta / denominator;
00328   covEtaPhi_ = numeratorEtaPhi / denominator;
00329   covPhiPhi_ = numeratorPhiPhi / denominator;
00330 }

void PFClusterShapeAlgo::fill5x5Map (  )  [private]

Definition at line 155 of file PFClusterShapeAlgo.cc.

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

Referenced by makeClusterShape().

00156 {
00157   // first get a navigator to the central element
00158   CaloNavigator<DetId> position = CaloNavigator<DetId>(eMaxId_, topoVector[topoIndex]);
00159 
00160   meanPosition_ = math::XYZVector(0.0, 0.0, 0.0);
00161   totalE_ = 0;
00162 
00163   for (int i = 0; i < 5; ++i)
00164     {
00165       for (int j = 0; j < 5; ++j)
00166         {
00167           position.home();
00168           position.offsetBy(i - 2, j - 2);
00169 
00170           RecHitWithFraction newEntry;
00171           newEntry.detId = DetId(0);
00172           newEntry.energy = 0.0;
00173           newEntry.position = math::XYZVector(0, 0, 0);
00174 
00175           if (*position != DetId(0)) // if this is a valid detId...
00176             {
00177               // ...find the corresponding PFRecHit index
00178               const int index = findPFRHIndexFromDetId((*position).rawId());
00179 
00180               if (index >= 0) // if a PFRecHit exists for this detId
00181                 {
00182                   double fraction = 1.0;
00183                   if (useFractions_) // if the algorithm should use fractions
00184                     { 
00185                       fraction = 0.0;
00186                       const reco::PFRecHitFraction * fraction_p = getFractionFromDetId(*position);
00187                       if (fraction_p) { fraction = fraction_p->fraction(); }
00188                     }
00189 
00190                   const reco::PFRecHitRef rhRef(currentRecHit_v_p, index);
00191                   const math::XYZVector crystalPosition(rhRef->position());
00192                   const double energyFraction =  rhRef->energy() * fraction;
00193 
00194                   newEntry.detId = *position;
00195                   newEntry.energy = energyFraction;
00196                   newEntry.position = crystalPosition;
00197 
00198                   meanPosition_ = meanPosition_ + crystalPosition * energyFraction; 
00199                   totalE_ += energyFraction;
00200                 }
00201             }
00202           map5x5[i][j] = newEntry;
00203         }
00204     }
00205   meanPosition_ /= totalE_;
00206 }

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().

00243 {
00244   std::map<double, Direction> directionMap;
00245 
00246   directionMap[addMapEnergies(-1, +0, -1, +0)] = SE;
00247   directionMap[addMapEnergies(-1, +0, +0, +1)] = NE;
00248   directionMap[addMapEnergies(+0, +1, -1, +0)] = SW;
00249   directionMap[addMapEnergies(+0, +1, +0, +1)] = NW;
00250 
00251   const std::map<double, Direction>::reverse_iterator eMaxDir_it = directionMap.rbegin();
00252 
00253   eMaxDir = eMaxDir_it->second;
00254 
00255   e2x2_ = eMaxDir_it->first;
00256 }

void PFClusterShapeAlgo::find_e2x5Bottom (  )  [private]

Definition at line 232 of file PFClusterShapeAlgo.cc.

References addMapEnergies(), and e2x5Bottom_.

Referenced by makeClusterShape().

00232 { e2x5Bottom_ = addMapEnergies(+1, +2, -2, +2); }

void PFClusterShapeAlgo::find_e2x5Left (  )  [private]

Definition at line 230 of file PFClusterShapeAlgo.cc.

References addMapEnergies(), and e2x5Left_.

Referenced by makeClusterShape().

00230 { e2x5Left_   = addMapEnergies(-2, +2, -2, -1); }

void PFClusterShapeAlgo::find_e2x5Right (  )  [private]

Definition at line 229 of file PFClusterShapeAlgo.cc.

References addMapEnergies(), and e2x5Right_.

Referenced by makeClusterShape().

00229 { e2x5Right_  = addMapEnergies(-2, +2, +1, +2); }

void PFClusterShapeAlgo::find_e2x5Top (  )  [private]

Definition at line 231 of file PFClusterShapeAlgo.cc.

References addMapEnergies(), and e2x5Top_.

Referenced by makeClusterShape().

00231 { e2x5Top_    = addMapEnergies(-2, -1, -2, +2); }

void PFClusterShapeAlgo::find_e3x2 (  )  [private]

Definition at line 258 of file PFClusterShapeAlgo.cc.

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

Referenced by makeClusterShape().

00259 {
00260   // Find the direction of the highest energy neighbour
00261   std::map<double, Direction> directionMap;
00262   directionMap[map5x5[2][3].energy] = N; 
00263   directionMap[map5x5[2][1].energy] = S;
00264   directionMap[map5x5[1][2].energy] = E;
00265   directionMap[map5x5[3][2].energy] = W;
00266   // Maps are sorted in ascending order - get the last element
00267   const Direction dir = directionMap.rbegin()->second;
00268 
00269   if (dir == N) 
00270     {
00271       e3x2_ = addMapEnergies(-1, +1, +0, +1);
00272       const double numerator   = map5x5[3][2].energy + map5x5[1][2].energy;
00273       const double denominator = map5x5[1][3].energy + map5x5[3][3].energy + 0.5;
00274       e3x2Ratio_ = numerator / denominator;
00275     }
00276   else if (dir == S)
00277     {
00278       e3x2_ = addMapEnergies(-1, +1, -1, +0);
00279       const double numerator   = map5x5[3][2].energy + map5x5[1][2].energy;
00280       const double denominator = map5x5[1][1].energy + map5x5[3][1].energy + 0.5;
00281       e3x2Ratio_ = numerator / denominator;
00282     }
00283   else if (dir == W)
00284     {
00285       e3x2_ = addMapEnergies(+0, +1, -1, +1);
00286       const double numerator   = map5x5[2][3].energy + map5x5[2][1].energy;
00287       const double denominator = map5x5[3][3].energy + map5x5[3][1].energy + 0.5;
00288       e3x2Ratio_ = numerator / denominator;
00289     }
00290   else if (dir == E)
00291     {
00292       e3x2_ = addMapEnergies(-1, +0, -1, +1);
00293       const double numerator   = map5x5[2][3].energy + map5x5[2][1].energy;
00294       const double denominator = map5x5[1][1].energy + map5x5[1][3].energy + 0.5;
00295       e3x2Ratio_ = numerator / denominator;
00296     }
00297 }

void PFClusterShapeAlgo::find_e3x3 (  )  [private]

Definition at line 227 of file PFClusterShapeAlgo.cc.

References addMapEnergies(), and e3x3_.

Referenced by makeClusterShape().

00227 { e3x3_       = addMapEnergies(-1, +1, -1, +1); }

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().

00235 {
00236   if (eMaxDir == SE) { e4x4_ = addMapEnergies(-2, +1, -2, +1); return; }
00237   if (eMaxDir == NE) { e4x4_ = addMapEnergies(-2, +1, -1, +2); return; }
00238   if (eMaxDir == SW) { e4x4_ = addMapEnergies(-1, +2, -2, +1); return; }
00239   if (eMaxDir == NW) { e4x4_ = addMapEnergies(-1, +2, -1, +2); return; }
00240 }

void PFClusterShapeAlgo::find_e5x5 (  )  [private]

Definition at line 228 of file PFClusterShapeAlgo.cc.

References addMapEnergies(), and e5x5_.

Referenced by makeClusterShape().

00228 { e5x5_       = addMapEnergies(-2, +2, -2, +2); }

void PFClusterShapeAlgo::find_eMax_e2nd (  )  [private]

Definition at line 96 of file PFClusterShapeAlgo.cc.

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

Referenced by makeClusterShape().

00097 {
00098   std::map<double, DetId> energyMap;
00099 
00100   // First get the RecHitFractions:
00101   const std::vector<reco::PFRecHitFraction> & fraction_v = currentCluster_p->recHitFractions();
00102   // For every one of them...
00103   for (std::vector<reco::PFRecHitFraction>::const_iterator it = fraction_v.begin(); it != fraction_v.end(); ++it)
00104     {
00105       // ...find the corresponding rechit:
00106       // const reco::PFRecHit & rechit = (*currentRecHit_v_p)[it->recHitIndex()];
00107       const reco::PFRecHitRef rechit = it->recHitRef();
00108       // ...and DetId:
00109       const DetId rechitDetId = DetId(rechit->detId());
00110       // Make the new Pair and put it in the map:
00111       energyMap[rechit->energy()] = rechitDetId;
00112     }
00113   // maps are sorted in ascending order so get the last two elements:
00114   std::map<double, DetId>::reverse_iterator it = energyMap.rbegin();
00115   eMax_   = it->first;
00116   eMaxId_ = it->second;
00117   it++;
00118   e2nd_   = it->first;
00119   e2ndId_ = it->second;
00120 }

int PFClusterShapeAlgo::findPFRHIndexFromDetId ( unsigned int  id  )  [private]

Definition at line 122 of file PFClusterShapeAlgo.cc.

References currentRecHit_v_p, index, and k.

Referenced by fill5x5Map().

00123 {
00124   int index = -1; // need some negative number
00125   for (unsigned int k = 0; k < currentRecHit_v_p->size(); ++k)
00126     {
00127       if ((*currentRecHit_v_p)[k].detId() == id)
00128         {
00129           index = static_cast<int>(k);
00130           break;
00131         }
00132     }
00133   return index;
00134 }

const reco::PFRecHitFraction * PFClusterShapeAlgo::getFractionFromDetId ( const DetId id  )  [private]

Definition at line 137 of file PFClusterShapeAlgo.cc.

References currentCluster_p, and it.

Referenced by fill5x5Map().

00138 {
00139   const std::vector< reco::PFRecHitFraction > & fraction_v = currentCluster_p->recHitFractions();
00140   for (std::vector<reco::PFRecHitFraction>::const_iterator it = fraction_v.begin(); it != fraction_v.end(); ++it)
00141     {
00142       //const unsigned int rhIndex = it->recHitIndex();
00143       //reco::PFRecHitRef rh_p(currentRecHit_v_p, rhIndex);
00144       const reco::PFRecHitRef rh_p = it->recHitRef();
00145       const DetId rhDetId = DetId(rh_p->detId());
00146       if (rhDetId == id) 
00147         { 
00148           return &(*it); 
00149         }
00150     }
00151   return 0;
00152 }

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().

00061 {
00062   find_eMax_e2nd();
00063   fill5x5Map();
00064   
00065   find_e2x2();
00066   find_e3x2();
00067   find_e3x3();
00068   find_e4x4();
00069   find_e5x5();
00070   
00071   find_e2x5Right();
00072   find_e2x5Left();
00073   find_e2x5Top();
00074   find_e2x5Bottom();
00075   
00076   covariances();
00077   
00078   double dummyLAT = 0;
00079   double dummyEtaLAT = 0;
00080   double dummyPhiLAT = 0;
00081   double dummyA20 = 0;
00082   double dummyA42 = 0;
00083 
00084   std::vector<double> dummyEnergyBasketFractionEta_v;
00085   std::vector<double> dummyEnergyBasketFractionPhi_v;
00086 
00087   return reco::ClusterShape(covEtaEta_, covEtaPhi_, covPhiPhi_, 
00088                             eMax_, eMaxId_, e2nd_, e2ndId_,
00089                             e2x2_, e3x2_, e3x3_, e4x4_, e5x5_, 
00090                             e2x5Right_, e2x5Left_, e2x5Top_, e2x5Bottom_, e3x2Ratio_,
00091                             dummyLAT, dummyEtaLAT, dummyPhiLAT, dummyA20, dummyA42,
00092                             dummyEnergyBasketFractionEta_v, dummyEnergyBasketFractionPhi_v);
00093 }

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, i, makeClusterShape(), topoIndex, and topoVector.

Referenced by PFClusterShapeProducer::produce().

00020 {
00021   static const float etaEndOfBarrel = 1.497;
00022 
00023   topoVector.push_back(the_barrelTop_p);
00024   topoVector.push_back(the_endcapTop_p);
00025   geomVector.push_back(the_barrelGeo_p);
00026   geomVector.push_back(the_endcapGeo_p);
00027 
00028   reco::ClusterShapeCollection * shape_v_p = new reco::ClusterShapeCollection();
00029 
00030   currentRecHit_v_p = rechitHandle;
00031 
00032   for (unsigned int i = 0; i < clusterHandle->size(); ++i)
00033     {
00034       // Make each cluster the "current" cluster
00035       currentCluster_p = reco::PFClusterRef(clusterHandle, i);
00036       currentClusterIndex_ = i;
00037 
00038       // Find the right topology to use with this cluster
00039       topoIndex = BARREL;
00040       geomIndex = BARREL;
00041       const math::XYZVector currentClusterPos(currentCluster_p->position());
00042       if (fabs(currentClusterPos.eta()) > etaEndOfBarrel)
00043         {
00044           topoIndex = ENDCAP;
00045           geomIndex = ENDCAP;
00046         }
00047 
00048       // Create the clustershape and push it into the vector
00049       shape_v_p->push_back(makeClusterShape());
00050     }
00051 
00052   topoVector.clear();
00053   topoVector.clear();
00054   geomVector.clear();
00055   geomVector.clear();
00056 
00057   return shape_v_p;
00058 }


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().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:29:42 2009 for CMSSW by  doxygen 1.5.4