CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

PFRecHitProducerECAL Class Reference

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

#include <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, std::vector< reco::PFRecHit > &rechitsCleaned, edm::Event &, const edm::EventSetup &)
void ecalNeighbArray (const EcalBarrelGeometry &barrelGeom, const CaloSubdetectorTopology &barrelTopo, const EcalEndcapGeometry &endcapGeom, const CaloSubdetectorTopology &endcapTopo)
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
double threshCleaning_
bool timingCleaning_

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_, neighbourmapcalculated_, threshCleaning_, and timingCleaning_.

 : PFRecHitProducer(iConfig) {

  // access to the collections of rechits

  
  inputTagEcalRecHitsEB_ = 
    iConfig.getParameter<InputTag>("ecalRecHitsEB");

  inputTagEcalRecHitsEE_ = 
    iConfig.getParameter<InputTag>("ecalRecHitsEE");
  
  
  crossBarrelEndcapBorder_ =
    iConfig.getParameter<bool>("crossBarrelEndcapBorder");

  timingCleaning_ = 
    iConfig.getParameter<bool>("timing_Cleaning");

  threshCleaning_ = 
    iConfig.getParameter<double>("thresh_Cleaning");

  
  neighbourmapcalculated_ = false;
}
PFRecHitProducerECAL::~PFRecHitProducerECAL ( )

Definition at line 64 of file PFRecHitProducerECAL.cc.

{}

Member Function Documentation

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

Definition at line 252 of file PFRecHitProducerECAL.cc.

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

Referenced by PFRecHitProducerPS::createRecHits().

                                                                              {

  math::XYZVector position;
  math::XYZVector axis;

  const CaloCellGeometry *thisCell 
    = geom->getGeometry(detid);
  
  // find rechit geometry
  if(!thisCell) {
    LogError("PFRecHitProducerECAL")
      <<"warning detid "<<detid.rawId()
      <<" not found in geometry"<<endl;
    return 0;
  }
  
  position.SetCoordinates ( thisCell->getPosition().x(),
                            thisCell->getPosition().y(),
                            thisCell->getPosition().z() );

  
  
  // the axis vector is the difference 
  const TruncatedPyramid* pyr 
    = dynamic_cast< const TruncatedPyramid* > (thisCell);    
  if( pyr ) {
    axis.SetCoordinates( pyr->getPosition(1).x(), 
                         pyr->getPosition(1).y(), 
                         pyr->getPosition(1).z() ); 
    
    math::XYZVector axis0( pyr->getPosition(0).x(), 
                           pyr->getPosition(0).y(), 
                           pyr->getPosition(0).z() );
    
    axis -= axis0;    
  }
  else return 0;

//   if( !geomfound ) {
//     LogError("PFRecHitProducerECAL")<<"cannor find geometry for detid "
//                               <<detid.rawId()<<" in layer "<<layer<<endl;
//     return 0; // geometry not found, skip rechit
//   }
  
  
  reco::PFRecHit *rh 
    = new reco::PFRecHit( detid.rawId(), layer, 
                          energy, 
                          position.x(), position.y(), position.z(), 
                          axis.x(), axis.y(), axis.z() ); 


  const CaloCellGeometry::CornersVec& corners = thisCell->getCorners();
  assert( corners.size() == 8 );

  rh->setNECorner( corners[0].x(), corners[0].y(),  corners[0].z() );
  rh->setSECorner( corners[1].x(), corners[1].y(),  corners[1].z() );
  rh->setSWCorner( corners[2].x(), corners[2].y(),  corners[2].z() );
  rh->setNWCorner( corners[3].x(), corners[3].y(),  corners[3].z() );

  return rh;
}
void PFRecHitProducerECAL::createRecHits ( std::vector< reco::PFRecHit > &  rechits,
std::vector< reco::PFRecHit > &  rechitsCleaned,
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 422 of file PFRecHitProducerECAL.cc.

References EAST, DetId::Ecal, EcalBarrel, EcalEndcap, CaloSubdetectorGeometry::getValidDetIds(), CaloSubdetectorTopology::getWindow(), EEDetId::hashedIndex(), EBDetId::hashedIndex(), recoMuon::in, LogDebug, neighbourmapcalculated_, neighboursEB_, neighboursEE_, NORTH, NORTHEAST, NORTHWEST, findQualityFiles::size, SOUTH, SOUTHEAST, SOUTHWEST, ntuplemaker::status, stdmove(), and WEST.

Referenced by PFRecHitProducerPS::createRecHits().

                                                                         {
  

  static const CaloDirection orderedDir[8]={SOUTHWEST,
                                            SOUTH,
                                            SOUTHEAST,
                                            WEST,
                                            EAST,
                                            NORTHWEST,
                                            NORTH,
                                            NORTHEAST};

  const unsigned nbarrel = 62000;
  // Barrel first. The hashed index runs from 0 to 61199
  neighboursEB_.resize(nbarrel);
  
  //std::cout << " Building the array of neighbours (barrel) " ;

  const std::vector<DetId>& vec(barrelGeom.getValidDetIds(DetId::Ecal,
                                                          EcalBarrel));
  unsigned size=vec.size();    
  for(unsigned ic=0; ic<size; ++ic) 
    {
      // We get the 9 cells in a square. 
      std::vector<DetId> neighbours(barrelTopo.getWindow(vec[ic],3,3));
      //      std::cout << " Cell " << EBDetId(vec[ic]) << std::endl;
      unsigned nneighbours=neighbours.size();

      unsigned hashedindex=EBDetId(vec[ic]).hashedIndex();
      if(hashedindex>=nbarrel)
        {
          LogDebug("CaloGeometryTools")  << " Array overflow " << std::endl;
        }


      // If there are 9 cells, it is easy, and this order is know:
//      6  7  8
//      3  4  5 
//      0  1  2   (0 = SOUTHWEST)

      if(nneighbours==9)
        {
          neighboursEB_[hashedindex].reserve(8);
          for(unsigned in=0;in<nneighbours;++in)
            {
              // remove the centre
              if(neighbours[in]!=vec[ic]) 
                {
                  neighboursEB_[hashedindex].push_back(neighbours[in]);
                  //          std::cout << " Neighbour " << in << " " << EBDetId(neighbours[in]) << std::endl;
                }
            }
        }
      else
        {
          DetId central(vec[ic]);
          neighboursEB_[hashedindex].resize(8,DetId(0));
          for(unsigned idir=0;idir<8;++idir)
            {
              DetId testid=central;
              bool status=stdmove(testid,orderedDir[idir],
                                  barrelTopo, endcapTopo,
                                  barrelGeom, endcapGeom);
              if(status) neighboursEB_[hashedindex][idir]=testid;
            }

        }
    }

  // Moved to the endcap

  //  std::cout << " done " << size << std::endl;
  //  std::cout << " Building the array of neighbours (endcap) " ;

//  vec.clear();
  const std::vector<DetId>& vecee=endcapGeom.getValidDetIds(DetId::Ecal,EcalEndcap);
  size=vecee.size();    
  // There are some holes in the hashedIndex for the EE. Hence the array is bigger than the number
  // of crystals
  const unsigned nendcap=19960;

  neighboursEE_.resize(nendcap);
  for(unsigned ic=0; ic<size; ++ic) 
    {
      // We get the 9 cells in a square. 
      std::vector<DetId> neighbours(endcapTopo.getWindow(vecee[ic],3,3));
      unsigned nneighbours=neighbours.size();
      // remove the centre
      unsigned hashedindex=EEDetId(vecee[ic]).hashedIndex();
      
      if(hashedindex>=nendcap)
        {
          LogDebug("CaloGeometryTools")  << " Array overflow " << std::endl;
        }

      if(nneighbours==9)
        {
          neighboursEE_[hashedindex].reserve(8);
          for(unsigned in=0;in<nneighbours;++in)
            {     
              // remove the centre
              if(neighbours[in]!=vecee[ic]) 
                {
                  neighboursEE_[hashedindex].push_back(neighbours[in]);
                }
            }
        }
      else
        {
          DetId central(vecee[ic]);
          neighboursEE_[hashedindex].resize(8,DetId(0));
          for(unsigned idir=0;idir<8;++idir)
            {
              DetId testid=central;
              bool status=stdmove(testid,orderedDir[idir],
                                  barrelTopo, endcapTopo,
                                  barrelGeom, endcapGeom);

              if(status) neighboursEE_[hashedindex][idir]=testid;
            }

        }
    }
  //  std::cout << " done " << size <<std::endl;
  neighbourmapcalculated_ = true;
}
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 322 of file PFRecHitProducerECAL.cc.

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

                                                                {
  

  const CaloCellGeometry *thisCell 
    = geom->getGeometry(detid);
  
  // find rechit geometry
  if(!thisCell) {
    LogError("PFRecHitProducerECAL")
      <<"warning detid "<<detid.rawId()
      <<" not found in geometry"<<endl;
    return false;
  }
  
  position.SetCoordinates ( thisCell->getPosition().x(),
                            thisCell->getPosition().y(),
                            thisCell->getPosition().z() );

  
  
  // the axis vector is the difference 
  const TruncatedPyramid* pyr 
    = dynamic_cast< const TruncatedPyramid* > (thisCell);    
  if( pyr ) {
    axis.SetCoordinates( pyr->getPosition(1).x(), 
                         pyr->getPosition(1).y(), 
                         pyr->getPosition(1).z() ); 
    
    math::XYZVector axis0( pyr->getPosition(0).x(), 
                           pyr->getPosition(0).y(), 
                           pyr->getPosition(0).z() );
    
    axis -= axis0;

    
    return true;
  }
  else return false;
}
void PFRecHitProducerECAL::findRecHitNeighboursECAL ( reco::PFRecHit rh,
const std::map< unsigned, unsigned > &  sortedHits 
) [private]

find rechit neighbours, using the hashed index

Definition at line 369 of file PFRecHitProducerECAL.cc.

References reco::PFRecHit::add4Neighbour(), reco::PFRecHit::add8Neighbour(), reco::PFRecHit::detId(), EAST, i, NORTH, NORTHEAST, NORTHWEST, DetId::rawId(), SOUTH, SOUTHEAST, SOUTHWEST, and WEST.

Referenced by PFRecHitProducerPS::createRecHits().

                                              {
  
  DetId center( rh.detId() );


  DetId north = move( center, NORTH );
  DetId northeast = move( center, NORTHEAST );
  DetId northwest = move( center, NORTHWEST ); 
  DetId south = move( center, SOUTH );  
  DetId southeast = move( center, SOUTHEAST );  
  DetId southwest = move( center, SOUTHWEST );  
  DetId east  = move( center, EAST );  
  DetId west  = move( center, WEST );  
    
  IDH i = sortedHits.find( north.rawId() );
  if(i != sortedHits.end() ) 
    rh.add4Neighbour( i->second );
  
  i = sortedHits.find( northeast.rawId() );
  if(i != sortedHits.end() ) 
    rh.add8Neighbour( i->second );
  
  i = sortedHits.find( south.rawId() );
  if(i != sortedHits.end() ) 
    rh.add4Neighbour( i->second );
    
  i = sortedHits.find( southwest.rawId() );
  if(i != sortedHits.end() ) 
    rh.add8Neighbour( i->second );
    
  i = sortedHits.find( east.rawId() );
  if(i != sortedHits.end() ) 
    rh.add4Neighbour( i->second );
    
  i = sortedHits.find( southeast.rawId() );
  if(i != sortedHits.end() ) 
    rh.add8Neighbour( i->second );
    
  i = sortedHits.find( west.rawId() );
  if(i != sortedHits.end() ) 
     rh.add4Neighbour( i->second );
   
  i = sortedHits.find( northwest.rawId() );
  if(i != sortedHits.end() ) 
    rh.add8Neighbour( i->second );
}
DetId PFRecHitProducerECAL::move ( DetId  cell,
const CaloDirection dir 
) const [private]

Definition at line 727 of file PFRecHitProducerECAL.cc.

References dir, EcalBarrel, EEDetId::hashedIndex(), EBDetId::hashedIndex(), neighbourmapcalculated_, neighboursEB_, neighboursEE_, NONE, query::result, and DetId::subdetId().

{  
  DetId originalcell = cell; 
  if(dir==NONE || cell==DetId(0)) return false;

  // Conversion CaloDirection and index in the table
  // CaloDirection :NONE,SOUTH,SOUTHEAST,SOUTHWEST,EAST,WEST, NORTHEAST,NORTHWEST,NORTH
  // Table : SOUTHWEST,SOUTH,SOUTHEAST,WEST,EAST,NORTHWEST,NORTH, NORTHEAST
  static const int calodirections[9]={-1,1,2,0,4,3,7,5,6};
    
  assert(neighbourmapcalculated_);

  DetId result = (originalcell.subdetId()==EcalBarrel) ? 
    neighboursEB_[EBDetId(originalcell).hashedIndex()][calodirections[dir]]:
    neighboursEE_[EEDetId(originalcell).hashedIndex()][calodirections[dir]];
  return result; 
}
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 634 of file PFRecHitProducerECAL.cc.

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

Referenced by ecalNeighbArray().

        {


  bool result; 

  if(dir==NORTH) {
    result = stdsimplemove(cell,NORTH, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
    return result;
  }
  else if(dir==SOUTH) {
    result = stdsimplemove(cell,SOUTH, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
    return result;
  }
  else if(dir==EAST) {
    result = stdsimplemove(cell,EAST, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
    return result;
  }
  else if(dir==WEST) {
    result = stdsimplemove(cell,WEST, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
    return result;
  }


  // One has to try both paths
  else if(dir==NORTHEAST)
    {
      result = stdsimplemove(cell,NORTH, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
      if(result)
        return stdsimplemove(cell,EAST, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
      else
        {
          result = stdsimplemove(cell,EAST, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
          if(result)
            return stdsimplemove(cell,NORTH, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
          else
            return false; 
        }
    }
  else if(dir==NORTHWEST)
    {
      result = stdsimplemove(cell,NORTH, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
      if(result)
        return stdsimplemove(cell,WEST, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
      else
        {
          result = stdsimplemove(cell,WEST, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
          if(result)
            return stdsimplemove(cell,NORTH, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
          else
            return false; 
        }
    }
  else if(dir == SOUTHEAST)
    {
      result = stdsimplemove(cell,SOUTH, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
      if(result)
        return stdsimplemove(cell,EAST, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
      else
        {
          result = stdsimplemove(cell,EAST, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
          if(result)
            return stdsimplemove(cell,SOUTH, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
          else
            return false; 
        }
    }
  else if(dir == SOUTHWEST)
    {
      result = stdsimplemove(cell,SOUTH, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
      if(result)
        return stdsimplemove(cell,WEST, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
      else
        {
          result = stdsimplemove(cell,SOUTH, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
          if(result)
            return stdsimplemove(cell,WEST, barrelTopo, endcapTopo, barrelGeom, endcapGeom );
          else
            return false; 
        }
    }
  cell = DetId(0);
  return false;
}
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 556 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().

        {

  std::vector<DetId> neighbours;

  // BARREL CASE 
  if(cell.subdetId()==EcalBarrel) {
    EBDetId ebDetId = cell;

    neighbours = barrelTopo.getNeighbours(ebDetId,dir);

    // first try to move according to the standard navigation
    if(neighbours.size()>0 && !neighbours[0].null()) {
      cell = neighbours[0];
      return true;
    }

    // failed.

    if(crossBarrelEndcapBorder_) {
      // are we on the outer ring ?
      const int ietaAbs ( ebDetId.ietaAbs() ) ; // abs value of ieta
      if( EBDetId::MAX_IETA == ietaAbs ) {
        // get ee nbrs for for end of barrel crystals  
        
        // yes we are
        const EcalBarrelGeometry::OrderedListOfEEDetId& 
          ol( * barrelGeom.getClosestEndcapCells( ebDetId ) ) ;
        
        // take closest neighbour on the other side, that is in the barrel.
        cell = *(ol.begin() );
        return true;
      }   
    }
  }

  // ENDCAP CASE 
  else if(cell.subdetId()==EcalEndcap) {

    EEDetId eeDetId = cell;

    neighbours= endcapTopo.getNeighbours(eeDetId,dir);

    if(neighbours.size()>0 && !neighbours[0].null()) {
      cell = neighbours[0];
      return true;
    }

    // failed.

    if(crossBarrelEndcapBorder_) {
      // are we on the outer ring ?
      const int iphi ( eeDetId.iPhiOuterRing() ) ;    
      if( iphi!= 0) {
        // yes we are
        const EcalEndcapGeometry::OrderedListOfEBDetId& 
          ol( * endcapGeom.getClosestBarrelCells( eeDetId ) ) ;
        
        // take closest neighbour on the other side, that is in the barrel.
        cell = *(ol.begin() );
        return true;
      }   
    }
  } 

  // everything failed 
  cell = DetId(0);
  return false;
}

Member Data Documentation

if true, navigation will cross the barrel-endcap border

Definition at line 117 of file PFRecHitProducerECAL.h.

Referenced by PFRecHitProducerECAL(), and stdsimplemove().

set to true in ecalNeighbArray

Definition at line 114 of file PFRecHitProducerECAL.h.

Referenced by PFRecHitProducerPS::createRecHits(), 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 108 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 111 of file PFRecHitProducerECAL.h.

Referenced by ecalNeighbArray(), and move().