CMS 3D CMS Logo

PFRecHitProducerECAL Class Reference

Producer for particle flow rechits (PFRecHit) in ECAL. More...

#include <RecoParticleFlow/PFClusterProducer/interface/PFRecHitProducerECAL.h>

Inheritance diagram for PFRecHitProducerECAL:

PFRecHitProducer edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Public Member Functions

 PFRecHitProducerECAL (const edm::ParameterSet &)
 ~PFRecHitProducerECAL ()

Private Member Functions

reco::PFRecHitcreateEcalRecHit (const DetId &detid, double energy, PFLayer::Layer layer, const CaloSubdetectorGeometry *geom)
void createRecHits (std::vector< reco::PFRecHit > &rechits, edm::Event &, const edm::EventSetup &)
 gets ecal barrel and endcap rechits, translate them to PFRecHits, which are stored in the rechits vector
void ecalNeighbArray (const EcalBarrelGeometry &barrelGeom, const CaloSubdetectorTopology &barrelTopo, const EcalEndcapGeometry &endcapGeom, const CaloSubdetectorTopology &endcapTopo)
 fill the vectors neighboursEB_ and neighboursEE_ which keep track of the neighbours of each rechit.
bool findEcalRecHitGeometry (const DetId &detid, const CaloSubdetectorGeometry *geom, math::XYZVector &position, math::XYZVector &axis)
 find the position and the axis of the cell for a given rechit
void findRecHitNeighboursECAL (reco::PFRecHit &rh, const std::map< unsigned, unsigned > &sortedHits)
 find rechit neighbours, using the hashed index
DetId move (DetId cell, const CaloDirection &dir) const
bool stdmove (DetId &cell, const CaloDirection &dir, const CaloSubdetectorTopology &barrelTopo, const CaloSubdetectorTopology &endcapTopo, const EcalBarrelGeometry &barrelGeom, const EcalEndcapGeometry &endcapGeom) const
bool stdsimplemove (DetId &cell, const CaloDirection &dir, const CaloSubdetectorTopology &barrelTopo, const CaloSubdetectorTopology &endcapTopo, const EcalBarrelGeometry &barrelGeom, const EcalEndcapGeometry &endcapGeom) const

Private Attributes

bool crossBarrelEndcapBorder_
 if true, navigation will cross the barrel-endcap border
edm::InputTag inputTagEcalRecHitsEB_
edm::InputTag inputTagEcalRecHitsEE_
bool neighbourmapcalculated_
 set to true in ecalNeighbArray
std::vector< std::vector< DetId > > neighboursEB_
 for each ecal barrel rechit, keep track of the neighbours
std::vector< std::vector< DetId > > neighboursEE_
 for each ecal endcap rechit, keep track of the neighbours


Detailed Description

Producer for particle flow rechits (PFRecHit) in ECAL.

Author:
Colin Bernet
Date:
february 2008

Definition at line 41 of file PFRecHitProducerECAL.h.


Constructor & Destructor Documentation

PFRecHitProducerECAL::PFRecHitProducerECAL ( const edm::ParameterSet iConfig  )  [explicit]

Definition at line 36 of file PFRecHitProducerECAL.cc.

References crossBarrelEndcapBorder_, edm::ParameterSet::getParameter(), inputTagEcalRecHitsEB_, inputTagEcalRecHitsEE_, and neighbourmapcalculated_.

00037  : PFRecHitProducer(iConfig) {
00038 
00039   // access to the collections of rechits
00040 
00041   
00042   inputTagEcalRecHitsEB_ = 
00043     iConfig.getParameter<InputTag>("ecalRecHitsEB");
00044 
00045   inputTagEcalRecHitsEE_ = 
00046     iConfig.getParameter<InputTag>("ecalRecHitsEE");
00047   
00048   
00049   crossBarrelEndcapBorder_ =
00050     iConfig.getParameter<bool>("crossBarrelEndcapBorder");
00051 
00052   
00053   neighbourmapcalculated_ = false;
00054 }

PFRecHitProducerECAL::~PFRecHitProducerECAL (  ) 

Definition at line 58 of file PFRecHitProducerECAL.cc.

00058 {}


Member Function Documentation

reco::PFRecHit * PFRecHitProducerECAL::createEcalRecHit ( const DetId detid,
double  energy,
PFLayer::Layer  layer,
const CaloSubdetectorGeometry geom 
) [private]

Definition at line 215 of file PFRecHitProducerECAL.cc.

References lat::endl(), CaloCellGeometry::getCorners(), CaloSubdetectorGeometry::getGeometry(), CaloCellGeometry::getPosition(), TruncatedPyramid::getPosition(), pyr, DetId::rawId(), reco::PFRecHit::setNECorner(), reco::PFRecHit::setNWCorner(), reco::PFRecHit::setSECorner(), reco::PFRecHit::setSWCorner(), PV3DBase< T, PVType, FrameType >::x(), x, PV3DBase< T, PVType, FrameType >::y(), y, PV3DBase< T, PVType, FrameType >::z(), and z.

00218                                                                               {
00219 
00220   math::XYZVector position;
00221   math::XYZVector axis;
00222 
00223   const CaloCellGeometry *thisCell 
00224     = geom->getGeometry(detid);
00225   
00226   // find rechit geometry
00227   if(!thisCell) {
00228     LogError("PFRecHitProducerECAL")
00229       <<"warning detid "<<detid.rawId()
00230       <<" not found in geometry"<<endl;
00231     return 0;
00232   }
00233   
00234   position.SetCoordinates ( thisCell->getPosition().x(),
00235                             thisCell->getPosition().y(),
00236                             thisCell->getPosition().z() );
00237 
00238   
00239   
00240   // the axis vector is the difference 
00241   const TruncatedPyramid* pyr 
00242     = dynamic_cast< const TruncatedPyramid* > (thisCell);    
00243   if( pyr ) {
00244     axis.SetCoordinates( pyr->getPosition(1).x(), 
00245                          pyr->getPosition(1).y(), 
00246                          pyr->getPosition(1).z() ); 
00247     
00248     math::XYZVector axis0( pyr->getPosition(0).x(), 
00249                            pyr->getPosition(0).y(), 
00250                            pyr->getPosition(0).z() );
00251     
00252     axis -= axis0;    
00253   }
00254   else return 0;
00255 
00256 //   if( !geomfound ) {
00257 //     LogError("PFRecHitProducerECAL")<<"cannor find geometry for detid "
00258 //                               <<detid.rawId()<<" in layer "<<layer<<endl;
00259 //     return 0; // geometry not found, skip rechit
00260 //   }
00261   
00262   
00263   reco::PFRecHit *rh 
00264     = new reco::PFRecHit( detid.rawId(), layer, 
00265                           energy, 
00266                           position.x(), position.y(), position.z(), 
00267                           axis.x(), axis.y(), axis.z() ); 
00268 
00269 
00270   const CaloCellGeometry::CornersVec& corners = thisCell->getCorners();
00271   assert( corners.size() == 8 );
00272 
00273   rh->setNECorner( corners[0].x(), corners[0].y(),  corners[0].z() );
00274   rh->setSECorner( corners[1].x(), corners[1].y(),  corners[1].z() );
00275   rh->setSWCorner( corners[2].x(), corners[2].y(),  corners[2].z() );
00276   rh->setNWCorner( corners[3].x(), corners[3].y(),  corners[3].z() );
00277 
00278   return rh;
00279 }

void PFRecHitProducerECAL::createRecHits ( std::vector< reco::PFRecHit > &  rechits,
edm::Event ,
const edm::EventSetup  
) [private, virtual]

gets ecal barrel and endcap rechits, translate them to PFRecHits, which are stored in the rechits vector

Implements PFRecHitProducer.

void PFRecHitProducerECAL::ecalNeighbArray ( const EcalBarrelGeometry barrelGeom,
const CaloSubdetectorTopology barrelTopo,
const EcalEndcapGeometry endcapGeom,
const CaloSubdetectorTopology endcapTopo 
) [private]

fill the vectors neighboursEB_ and neighboursEE_ which keep track of the neighbours of each rechit.

to be called at the beginning of the run

Definition at line 385 of file PFRecHitProducerECAL.cc.

References EAST, DetId::Ecal, EcalBarrel, EcalEndcap, lat::endl(), CaloSubdetectorGeometry::getValidDetIds(), CaloSubdetectorTopology::getWindow(), in, LogDebug, neighbourmapcalculated_, neighboursEB_, neighboursEE_, NORTH, NORTHEAST, NORTHWEST, size, SOUTH, SOUTHEAST, SOUTHWEST, StDecayID::status, stdmove(), and WEST.

00389                                                                          {
00390   
00391 
00392   static const CaloDirection orderedDir[8]={SOUTHWEST,
00393                                             SOUTH,
00394                                             SOUTHEAST,
00395                                             WEST,
00396                                             EAST,
00397                                             NORTHWEST,
00398                                             NORTH,
00399                                             NORTHEAST};
00400 
00401   const unsigned nbarrel = 62000;
00402   // Barrel first. The hashed index runs from 0 to 61199
00403   neighboursEB_.resize(nbarrel);
00404   
00405   //std::cout << " Building the array of neighbours (barrel) " ;
00406 
00407   std::vector<DetId> vec(barrelGeom.getValidDetIds(DetId::Ecal,
00408                                                    EcalBarrel));
00409   unsigned size=vec.size();    
00410   for(unsigned ic=0; ic<size; ++ic) 
00411     {
00412       // We get the 9 cells in a square. 
00413       std::vector<DetId> neighbours(barrelTopo.getWindow(vec[ic],3,3));
00414       //      std::cout << " Cell " << EBDetId(vec[ic]) << std::endl;
00415       unsigned nneighbours=neighbours.size();
00416 
00417       unsigned hashedindex=EBDetId(vec[ic]).hashedIndex();
00418       if(hashedindex>=nbarrel)
00419         {
00420           LogDebug("CaloGeometryTools")  << " Array overflow " << std::endl;
00421         }
00422 
00423 
00424       // If there are 9 cells, it is easy, and this order is know:
00425 //      6  7  8
00426 //      3  4  5 
00427 //      0  1  2   (0 = SOUTHWEST)
00428 
00429       if(nneighbours==9)
00430         {
00431           neighboursEB_[hashedindex].reserve(8);
00432           for(unsigned in=0;in<nneighbours;++in)
00433             {
00434               // remove the centre
00435               if(neighbours[in]!=vec[ic]) 
00436                 {
00437                   neighboursEB_[hashedindex].push_back(neighbours[in]);
00438                   //          std::cout << " Neighbour " << in << " " << EBDetId(neighbours[in]) << std::endl;
00439                 }
00440             }
00441         }
00442       else
00443         {
00444           DetId central(vec[ic]);
00445           neighboursEB_[hashedindex].resize(8,DetId(0));
00446           for(unsigned idir=0;idir<8;++idir)
00447             {
00448               DetId testid=central;
00449               bool status=stdmove(testid,orderedDir[idir],
00450                                   barrelTopo, endcapTopo,
00451                                   barrelGeom, endcapGeom);
00452               if(status) neighboursEB_[hashedindex][idir]=testid;
00453             }
00454 
00455         }
00456     }
00457 
00458   // Moved to the endcap
00459 
00460   //  std::cout << " done " << size << std::endl;
00461   //  std::cout << " Building the array of neighbours (endcap) " ;
00462 
00463   vec.clear();
00464   vec=endcapGeom.getValidDetIds(DetId::Ecal,EcalEndcap);
00465   size=vec.size();    
00466   // There are some holes in the hashedIndex for the EE. Hence the array is bigger than the number
00467   // of crystals
00468   const unsigned nendcap=19960;
00469 
00470   neighboursEE_.resize(nendcap);
00471   for(unsigned ic=0; ic<size; ++ic) 
00472     {
00473       // We get the 9 cells in a square. 
00474       std::vector<DetId> neighbours(endcapTopo.getWindow(vec[ic],3,3));
00475       unsigned nneighbours=neighbours.size();
00476       // remove the centre
00477       unsigned hashedindex=EEDetId(vec[ic]).hashedIndex();
00478       
00479       if(hashedindex>=nendcap)
00480         {
00481           LogDebug("CaloGeometryTools")  << " Array overflow " << std::endl;
00482         }
00483 
00484       if(nneighbours==9)
00485         {
00486           neighboursEE_[hashedindex].reserve(8);
00487           for(unsigned in=0;in<nneighbours;++in)
00488             {     
00489               // remove the centre
00490               if(neighbours[in]!=vec[ic]) 
00491                 {
00492                   neighboursEE_[hashedindex].push_back(neighbours[in]);
00493                 }
00494             }
00495         }
00496       else
00497         {
00498           DetId central(vec[ic]);
00499           neighboursEE_[hashedindex].resize(8,DetId(0));
00500           for(unsigned idir=0;idir<8;++idir)
00501             {
00502               DetId testid=central;
00503               bool status=stdmove(testid,orderedDir[idir],
00504                                   barrelTopo, endcapTopo,
00505                                   barrelGeom, endcapGeom);
00506 
00507               if(status) neighboursEE_[hashedindex][idir]=testid;
00508             }
00509 
00510         }
00511     }
00512   //  std::cout << " done " << size <<std::endl;
00513   neighbourmapcalculated_ = true;
00514 }

bool PFRecHitProducerECAL::findEcalRecHitGeometry ( const DetId detid,
const CaloSubdetectorGeometry geom,
math::XYZVector position,
math::XYZVector axis 
) [private]

find the position and the axis of the cell for a given rechit

Definition at line 285 of file PFRecHitProducerECAL.cc.

References lat::endl(), CaloSubdetectorGeometry::getGeometry(), CaloCellGeometry::getPosition(), TruncatedPyramid::getPosition(), pyr, DetId::rawId(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

00288                                                                 {
00289   
00290 
00291   const CaloCellGeometry *thisCell 
00292     = geom->getGeometry(detid);
00293   
00294   // find rechit geometry
00295   if(!thisCell) {
00296     LogError("PFRecHitProducerECAL")
00297       <<"warning detid "<<detid.rawId()
00298       <<" not found in geometry"<<endl;
00299     return false;
00300   }
00301   
00302   position.SetCoordinates ( thisCell->getPosition().x(),
00303                             thisCell->getPosition().y(),
00304                             thisCell->getPosition().z() );
00305 
00306   
00307   
00308   // the axis vector is the difference 
00309   const TruncatedPyramid* pyr 
00310     = dynamic_cast< const TruncatedPyramid* > (thisCell);    
00311   if( pyr ) {
00312     axis.SetCoordinates( pyr->getPosition(1).x(), 
00313                          pyr->getPosition(1).y(), 
00314                          pyr->getPosition(1).z() ); 
00315     
00316     math::XYZVector axis0( pyr->getPosition(0).x(), 
00317                            pyr->getPosition(0).y(), 
00318                            pyr->getPosition(0).z() );
00319     
00320     axis -= axis0;
00321 
00322     
00323     return true;
00324   }
00325   else return false;
00326 }

void PFRecHitProducerECAL::findRecHitNeighboursECAL ( reco::PFRecHit rh,
const std::map< unsigned, unsigned > &  sortedHits 
) [private]

find rechit neighbours, using the hashed index

DetId PFRecHitProducerECAL::move ( DetId  cell,
const CaloDirection dir 
) const [private]

Definition at line 690 of file PFRecHitProducerECAL.cc.

References EcalBarrel, neighbourmapcalculated_, neighboursEB_, neighboursEE_, NONE, HLT_VtxMuL3::result, and DetId::subdetId().

00692 {  
00693   DetId originalcell = cell; 
00694   if(dir==NONE || cell==DetId(0)) return false;
00695 
00696   // Conversion CaloDirection and index in the table
00697   // CaloDirection :NONE,SOUTH,SOUTHEAST,SOUTHWEST,EAST,WEST, NORTHEAST,NORTHWEST,NORTH
00698   // Table : SOUTHWEST,SOUTH,SOUTHEAST,WEST,EAST,NORTHWEST,NORTH, NORTHEAST
00699   static const int calodirections[9]={-1,1,2,0,4,3,7,5,6};
00700     
00701   assert(neighbourmapcalculated_);
00702 
00703   DetId result = (originalcell.subdetId()==EcalBarrel) ? 
00704     neighboursEB_[EBDetId(originalcell).hashedIndex()][calodirections[dir]]:
00705     neighboursEE_[EEDetId(originalcell).hashedIndex()][calodirections[dir]];
00706   return result; 
00707 }

bool PFRecHitProducerECAL::stdmove ( DetId cell,
const CaloDirection dir,
const CaloSubdetectorTopology barrelTopo,
const CaloSubdetectorTopology endcapTopo,
const EcalBarrelGeometry barrelGeom,
const EcalEndcapGeometry endcapGeom 
) const [private]

Definition at line 597 of file PFRecHitProducerECAL.cc.

References EAST, NORTH, NORTHEAST, NORTHWEST, HLT_VtxMuL3::result, SOUTH, SOUTHEAST, SOUTHWEST, stdsimplemove(), and WEST.

Referenced by ecalNeighbArray().

00604         {
00605 
00606 
00607   bool result; 
00608 
00609   if(dir==NORTH) {
00610     result = stdsimplemove(cell,NORTH, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
00611     return result;
00612   }
00613   else if(dir==SOUTH) {
00614     result = stdsimplemove(cell,SOUTH, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
00615     return result;
00616   }
00617   else if(dir==EAST) {
00618     result = stdsimplemove(cell,EAST, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
00619     return result;
00620   }
00621   else if(dir==WEST) {
00622     result = stdsimplemove(cell,WEST, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
00623     return result;
00624   }
00625 
00626 
00627   // One has to try both paths
00628   else if(dir==NORTHEAST)
00629     {
00630       result = stdsimplemove(cell,NORTH, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
00631       if(result)
00632         return stdsimplemove(cell,EAST, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
00633       else
00634         {
00635           result = stdsimplemove(cell,EAST, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
00636           if(result)
00637             return stdsimplemove(cell,NORTH, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
00638           else
00639             return false; 
00640         }
00641     }
00642   else if(dir==NORTHWEST)
00643     {
00644       result = stdsimplemove(cell,NORTH, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
00645       if(result)
00646         return stdsimplemove(cell,WEST, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
00647       else
00648         {
00649           result = stdsimplemove(cell,WEST, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
00650           if(result)
00651             return stdsimplemove(cell,NORTH, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
00652           else
00653             return false; 
00654         }
00655     }
00656   else if(dir == SOUTHEAST)
00657     {
00658       result = stdsimplemove(cell,SOUTH, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
00659       if(result)
00660         return stdsimplemove(cell,EAST, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
00661       else
00662         {
00663           result = stdsimplemove(cell,EAST, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
00664           if(result)
00665             return stdsimplemove(cell,SOUTH, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
00666           else
00667             return false; 
00668         }
00669     }
00670   else if(dir == SOUTHWEST)
00671     {
00672       result = stdsimplemove(cell,SOUTH, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
00673       if(result)
00674         return stdsimplemove(cell,WEST, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
00675       else
00676         {
00677           result = stdsimplemove(cell,SOUTH, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
00678           if(result)
00679             return stdsimplemove(cell,WEST, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
00680           else
00681             return false; 
00682         }
00683     }
00684   cell = DetId(0);
00685   return false;
00686 }

bool PFRecHitProducerECAL::stdsimplemove ( DetId cell,
const CaloDirection dir,
const CaloSubdetectorTopology barrelTopo,
const CaloSubdetectorTopology endcapTopo,
const EcalBarrelGeometry barrelGeom,
const EcalEndcapGeometry endcapGeom 
) const [private]

Definition at line 519 of file PFRecHitProducerECAL.cc.

References crossBarrelEndcapBorder_, EcalBarrel, EcalEndcap, EcalEndcapGeometry::getClosestBarrelCells(), EcalBarrelGeometry::getClosestEndcapCells(), CaloSubdetectorTopology::getNeighbours(), EBDetId::ietaAbs(), EEDetId::iPhiOuterRing(), EBDetId::MAX_IETA, and DetId::subdetId().

Referenced by stdmove().

00525         {
00526 
00527   std::vector<DetId> neighbours;
00528 
00529   // BARREL CASE 
00530   if(cell.subdetId()==EcalBarrel) {
00531     EBDetId ebDetId = cell;
00532 
00533     neighbours = barrelTopo.getNeighbours(ebDetId,dir);
00534 
00535     // first try to move according to the standard navigation
00536     if(neighbours.size()>0 && !neighbours[0].null()) {
00537       cell = neighbours[0];
00538       return true;
00539     }
00540 
00541     // failed.
00542 
00543     if(crossBarrelEndcapBorder_) {
00544       // are we on the outer ring ?
00545       const int ietaAbs ( ebDetId.ietaAbs() ) ; // abs value of ieta
00546       if( EBDetId::MAX_IETA == ietaAbs ) {
00547         // get ee nbrs for for end of barrel crystals  
00548         
00549         // yes we are
00550         const EcalBarrelGeometry::OrderedListOfEEDetId& 
00551           ol( * barrelGeom.getClosestEndcapCells( ebDetId ) ) ;
00552         
00553         // take closest neighbour on the other side, that is in the barrel.
00554         cell = *(ol.begin() );
00555         return true;
00556       }   
00557     }
00558   }
00559 
00560   // ENDCAP CASE 
00561   else if(cell.subdetId()==EcalEndcap) {
00562 
00563     EEDetId eeDetId = cell;
00564 
00565     neighbours= endcapTopo.getNeighbours(eeDetId,dir);
00566 
00567     if(neighbours.size()>0 && !neighbours[0].null()) {
00568       cell = neighbours[0];
00569       return true;
00570     }
00571 
00572     // failed.
00573 
00574     if(crossBarrelEndcapBorder_) {
00575       // are we on the outer ring ?
00576       const int iphi ( eeDetId.iPhiOuterRing() ) ;    
00577       if( iphi!= 0) {
00578         // yes we are
00579         const EcalEndcapGeometry::OrderedListOfEBDetId& 
00580           ol( * endcapGeom.getClosestBarrelCells( eeDetId ) ) ;
00581         
00582         // take closest neighbour on the other side, that is in the barrel.
00583         cell = *(ol.begin() );
00584         return true;
00585       }   
00586     }
00587   } 
00588 
00589   // everything failed 
00590   cell = DetId(0);
00591   return false;
00592 }


Member Data Documentation

bool PFRecHitProducerECAL::crossBarrelEndcapBorder_ [private]

if true, navigation will cross the barrel-endcap border

Definition at line 116 of file PFRecHitProducerECAL.h.

Referenced by PFRecHitProducerECAL(), and stdsimplemove().

edm::InputTag PFRecHitProducerECAL::inputTagEcalRecHitsEB_ [private]

Definition at line 119 of file PFRecHitProducerECAL.h.

Referenced by PFRecHitProducerECAL().

edm::InputTag PFRecHitProducerECAL::inputTagEcalRecHitsEE_ [private]

Definition at line 120 of file PFRecHitProducerECAL.h.

Referenced by PFRecHitProducerECAL().

bool PFRecHitProducerECAL::neighbourmapcalculated_ [private]

set to true in ecalNeighbArray

Definition at line 113 of file PFRecHitProducerECAL.h.

Referenced by ecalNeighbArray(), move(), and PFRecHitProducerECAL().

std::vector<std::vector<DetId> > PFRecHitProducerECAL::neighboursEB_ [private]

for each ecal barrel rechit, keep track of the neighbours

Definition at line 107 of file PFRecHitProducerECAL.h.

Referenced by ecalNeighbArray(), and move().

std::vector<std::vector<DetId> > PFRecHitProducerECAL::neighboursEE_ [private]

for each ecal endcap rechit, keep track of the neighbours

Definition at line 110 of file PFRecHitProducerECAL.h.

Referenced by ecalNeighbArray(), and move().


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