CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2_patch1/src/RecoParticleFlow/PFClusterProducer/plugins/PFRecHitProducerHO.cc

Go to the documentation of this file.
00001 #include "RecoParticleFlow/PFClusterProducer/plugins/PFRecHitProducerHO.h"
00002 
00003 #include <memory>
00004 
00005 #include "DataFormats/ParticleFlowReco/interface/PFRecHit.h"
00006 
00007 #include "DataFormats/HcalRecHit/interface/HORecHit.h"
00008 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
00009 #include "DataFormats/HcalDetId/interface/HcalSubdetector.h"
00010 
00011 #include "DataFormats/DetId/interface/DetId.h"
00012 #include "DataFormats/EcalDetId/interface/EEDetId.h"
00013 #include "DataFormats/EcalDetId/interface/EBDetId.h"
00014 
00015 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00016 #include "FWCore/Framework/interface/ESHandle.h"
00017 #include "FWCore/Framework/interface/EventSetup.h"
00018 
00019 #include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
00020 #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
00021 #include "Geometry/CaloGeometry/interface/CaloCellGeometry.h"
00022 #include "Geometry/CaloGeometry/interface/TruncatedPyramid.h"
00023 #include "Geometry/Records/interface/CaloGeometryRecord.h"
00024 
00025 #include "Geometry/CaloTopology/interface/HcalTopology.h"
00026 #include "DataFormats/CaloTowers/interface/CaloTowerCollection.h"
00027 
00028 #include "Geometry/CaloTopology/interface/CaloTowerTopology.h"
00029 #include "RecoCaloTools/Navigation/interface/CaloTowerNavigator.h"
00030 
00031 #include "RecoLocalCalo/HcalRecAlgos/interface/HcalCaloFlagLabels.h"
00032 #include "RecoCaloTools/Navigation/interface/CaloNavigator.h"
00033 
00034 // Necessary includes for identify severity of flagged problems in HO rechits
00035 #include "RecoLocalCalo/HcalRecAlgos/interface/HcalCaloFlagLabels.h"
00036 #include "RecoLocalCalo/HcalRecAlgos/interface/HcalSeverityLevelComputer.h"
00037 #include "RecoLocalCalo/HcalRecAlgos/interface/HcalSeverityLevelComputerRcd.h"
00038 #include "CondFormats/HcalObjects/interface/HcalChannelQuality.h"
00039 
00040 using namespace std;
00041 using namespace edm;
00042 
00043 PFRecHitProducerHO::PFRecHitProducerHO(const edm::ParameterSet& iConfig)
00044   : PFRecHitProducer(iConfig)
00045 {
00046   
00047   // access to the collections of rechits
00048   inputTagHORecHits_ = 
00049     iConfig.getParameter<InputTag>("recHitsHO");
00050   
00051   HOMaxAllowedSev_ = iConfig.getParameter<int>("HOMaxAllowedSev");
00052   neighbourmapcalculated_ = false;
00053 }
00054 
00055 
00056 
00057 PFRecHitProducerHO::~PFRecHitProducerHO() {}
00058 
00059 void 
00060 PFRecHitProducerHO::createRecHits(vector<reco::PFRecHit>& rechits,
00061                                   vector<reco::PFRecHit>& rechitsCleaned,
00062                                   edm::Event& iEvent, 
00063                                   const edm::EventSetup& iSetup ) {
00064 
00065 
00066 
00067   // this map is necessary to find the rechit neighbours efficiently
00068   //C but I should think about using Florian's hashed index to do this.
00069   //C in which case the map might not be necessary anymore
00070   // 
00071   // the key of this map is detId. 
00072   // the value is the index in the rechits vector
00073   map<unsigned, unsigned > idSortedRecHits;
00074   //   typedef map<unsigned, unsigned >::iterator IDH;
00075   
00076   edm::ESHandle<CaloGeometry> geoHandle;
00077   iSetup.get<CaloGeometryRecord>().get(geoHandle);
00078   
00079   // get the HO geometry
00080   const CaloSubdetectorGeometry *hcalBarrelGeometry = 
00081     geoHandle->getSubdetectorGeometry(DetId::Hcal, HcalOuter);
00082   
00083   // get the HO topology
00084   edm::ESHandle<HcalTopology> hcalBarrelTopology;
00085   iSetup.get<IdealGeometryRecord>().get(hcalBarrelTopology);
00086   
00087   if(!neighbourmapcalculated_)
00088     hoNeighbArray( *hcalBarrelGeometry,
00089                    *hcalBarrelTopology);
00090 
00091   // Get Hcal Severity Level Computer, so that the severity of each rechit flag/status may be determined
00092   edm::ESHandle<HcalSeverityLevelComputer> hcalSevLvlComputerHndl;
00093   iSetup.get<HcalSeverityLevelComputerRcd>().get(hcalSevLvlComputerHndl);
00094   const HcalSeverityLevelComputer* hcalSevLvlComputer = hcalSevLvlComputerHndl.product();
00095 
00096   
00097   // get the HO rechits
00098   
00099   edm::Handle<HORecHitCollection> rhcHandle;
00100   
00101   
00102   bool found = iEvent.getByLabel(inputTagHORecHits_, 
00103                                  rhcHandle);
00104   
00105   if(!found) {
00106     
00107     ostringstream err;
00108     err<<"could not find rechits "<<inputTagHORecHits_;
00109     LogError("PFRecHitProducerHO")<<err.str()<<endl;
00110     
00111     throw cms::Exception( "MissingProduct", err.str());
00112   }
00113   else {
00114     assert( rhcHandle.isValid() );
00115     
00116     // process HO rechits
00117     for(unsigned i=0; i<rhcHandle->size(); i++) {
00118       
00119       const HORecHit& erh = (*rhcHandle)[i];
00120       const HcalDetId& detid = (HcalDetId)erh.detid();
00121       
00122       double energy = erh.energy();
00123       // uint32_t flag = erh.recoFlag();
00124       double time = erh.time();
00125       
00126       HcalSubdetector esd=(HcalSubdetector)detid.subdetId();
00127       
00128       if (esd != 3) continue;
00129       
00130       int hoeta=detid.ieta();
00131       if ( (abs(hoeta)<=4 && energy < thresh_Barrel_) || 
00132            (abs(hoeta)> 4 && energy < thresh_Endcap_) ) continue;
00133       
00134 
00135       // Get Channel Quality information for the given detID
00136       const HcalChannelStatus* theStatus = theHcalChStatus->getValues(detid);
00137       unsigned theStatusValue = theStatus->getValue();
00138       // Now get severity of problems for the given detID, based on the rechit flag word and the channel quality status value
00139       int hitSeverity=hcalSevLvlComputer->getSeverityLevel(detid, erh.flags(),theStatusValue);
00140     
00141       // Skip hits whose problems are more severe than max accept level.  In the future, allow for cleaning of such hits?
00142       // Note:  As of April 2012, by default, all HO hits in rings +/-1, +/-2 should be identified as either "remove from calotowers" or "remove from rechit collections" in the channel quality database, and thus should be rejected by this conditional statement.
00143       if (hitSeverity>HOMaxAllowedSev_) 
00144         {
00145           //std::cout <<"Rejecting HO hit HO("<<hoeta<<", "<<detid.iphi()<<", "<<detid.depth()<<std::endl;
00146           continue;
00147         }  
00148 
00149 
00150       
00151       reco::PFRecHit *pfrh = createHORecHit(detid, energy,  
00152                                             PFLayer::HCAL_BARREL2, // HO,
00153                                             hcalBarrelGeometry);
00154       
00155       if( !pfrh ) continue; // problem with this rechit. skip it
00156       
00157       pfrh->setRescale(time);
00158       
00159       rechits.push_back( *pfrh );
00160       delete pfrh;
00161       idSortedRecHits.insert( make_pair(detid.rawId(), rechits.size()-1 ) ); 
00162     }      
00163   }
00164   
00165   // do navigation
00166   for(unsigned i=0; i<rechits.size(); i++ ) {
00167     findRecHitNeighboursHO( rechits[i], *hcalBarrelTopology, idSortedRecHits ); 
00168   }
00169   
00170 } 
00171 
00172 
00173 reco::PFRecHit* 
00174 PFRecHitProducerHO::createHORecHit( const DetId& detid,
00175                                     double energy,
00176                                     PFLayer::Layer layer,
00177                                     const CaloSubdetectorGeometry* geom ) {
00178   
00179   math::XYZVector position;
00180   math::XYZVector axis;
00181   
00182   const CaloCellGeometry *thisCell 
00183     = geom->getGeometry(detid);
00184   
00185   // find rechit geometry
00186   if(!thisCell) {
00187     LogError("PFRecHitProducerHO")
00188       <<"warning detid "<<detid.rawId()
00189       <<" not found in geometry"<<endl;
00190     return 0;
00191   }
00192   
00193   double sclel0l1r=0.946; //384.8/406.6
00194   
00195   if (abs(thisCell->getPosition().z())>130) {
00196     position.SetCoordinates ( sclel0l1r*thisCell->getPosition().x(),
00197                               sclel0l1r*thisCell->getPosition().y(),
00198                               sclel0l1r*thisCell->getPosition().z() );   
00199     
00200   } else {
00201     position.SetCoordinates ( thisCell->getPosition().x(),
00202                               thisCell->getPosition().y(),
00203                               thisCell->getPosition().z() );
00204   }
00205   
00206   // the axis vector is the difference 
00207   
00208   //   const TruncatedPyramid* pyr 
00209   //     = dynamic_cast< const TruncatedPyramid* > (thisCell);    
00210   //   if( pyr ) {
00211   //     axis.SetCoordinates( pyr->getPosition(1).x(), 
00212   //                     pyr->getPosition(1).y(), 
00213   //                     pyr->getPosition(1).z() ); 
00214   
00215   //     math::XYZVector axis0( pyr->getPosition(0).x(), 
00216   //                       pyr->getPosition(0).y(), 
00217   //                       pyr->getPosition(0).z() );
00218   
00219   //     axis -= axis0;    
00220   //   }
00221   //   else return 0;
00222   
00223   axis = math::XYZVector(0,0,0);
00224   
00225   //   if( !geomfound ) {
00226   //     LogError("PFRecHitProducerHO")<<"cannor find geometry for detid "
00227   //                             <<detid.rawId()<<" in layer "<<layer<<endl;
00228   //     return 0; // geometry not found, skip rechit
00229   //   }
00230   
00231   
00232   reco::PFRecHit *rh 
00233     = new reco::PFRecHit( detid.rawId(), layer, 
00234                           energy, 
00235                           position.x(), position.y(), position.z(), 
00236                           axis.x(), axis.y(), axis.z() ); 
00237   
00238   const CaloCellGeometry::CornersVec& corners = thisCell->getCorners();
00239   assert( corners.size() == 8 );
00240   
00241   if (abs(corners[0].z())>130.0) {
00242     rh->setNECorner( sclel0l1r*corners[0].x(), sclel0l1r*corners[0].y(),  sclel0l1r*corners[0].z() );
00243     rh->setSECorner( sclel0l1r*corners[1].x(), sclel0l1r*corners[1].y(),  sclel0l1r*corners[1].z() );
00244     rh->setSWCorner( sclel0l1r*corners[2].x(), sclel0l1r*corners[2].y(),  sclel0l1r*corners[2].z() );
00245     rh->setNWCorner( sclel0l1r*corners[3].x(), sclel0l1r*corners[3].y(),  sclel0l1r*corners[3].z() );
00246     
00247   } else {
00248     rh->setNECorner( corners[0].x(), corners[0].y(),  corners[0].z() );
00249     rh->setSECorner( corners[1].x(), corners[1].y(),  corners[1].z() );
00250     rh->setSWCorner( corners[2].x(), corners[2].y(),  corners[2].z() );
00251     rh->setNWCorner( corners[3].x(), corners[3].y(),  corners[3].z() );
00252   }
00253   
00254   return rh;
00255 }
00256 
00257 
00258 bool
00259 PFRecHitProducerHO::findHORecHitGeometry(const DetId& detid, 
00260                                          const CaloSubdetectorGeometry* geom,
00261                                          math::XYZVector& position, 
00262                                          math::XYZVector& axis ) {
00263   
00264   
00265   const CaloCellGeometry *thisCell 
00266     = geom->getGeometry(detid);
00267   
00268   // find rechit geometry
00269   if(!thisCell) {
00270     LogError("PFRecHitProducerHO")
00271       <<"warning detid "<<detid.rawId()
00272       <<" not found in geometry"<<endl;
00273     return false;
00274   }
00275   
00276   position.SetCoordinates ( thisCell->getPosition().x(),
00277                             thisCell->getPosition().y(),
00278                             thisCell->getPosition().z() );
00279   
00280   
00281   
00282   // the axis vector is the difference 
00283   //   const TruncatedPyramid* pyr 
00284   //     = dynamic_cast< const TruncatedPyramid* > (thisCell);    
00285   //   if( pyr ) {
00286   //     axis.SetCoordinates( pyr->getPosition(1).x(), 
00287   //                     pyr->getPosition(1).y(), 
00288   //                     pyr->getPosition(1).z() ); 
00289   
00290   //     math::XYZVector axis0( pyr->getPosition(0).x(), 
00291   //                       pyr->getPosition(0).y(), 
00292   //                       pyr->getPosition(0).z() );
00293   
00294   //     axis -= axis0;
00295   
00296   
00297   //     return true;
00298   //   }
00299   
00300   axis = math::XYZVector(0,0,0);
00301   return true;
00302   
00303   //  else return false;
00304 }
00305 
00306 
00307 
00308 void 
00309 PFRecHitProducerHO::findRecHitNeighboursHO
00310 ( reco::PFRecHit& rh, 
00311   const HcalTopology& topo, 
00312   const map<unsigned,unsigned >& sortedHits ) {
00313   
00314   DetId center( rh.detId() );
00315   
00316   
00317   DetId north = move( center, topo, NORTH );
00318   DetId northeast = move( center, topo, NORTHEAST );
00319   DetId northwest = move( center, topo, NORTHWEST ); 
00320   DetId south = move( center, topo, SOUTH );  
00321   DetId southeast = move( center, topo, SOUTHEAST );  
00322   DetId southwest = move( center, topo, SOUTHWEST );  
00323   DetId east  = move( center, topo, EAST );  
00324   DetId west  = move( center, topo, WEST );  
00325   
00326   IDH i = sortedHits.find( north.rawId() );
00327   if(i != sortedHits.end() ) 
00328     rh.add4Neighbour( i->second );
00329   
00330   i = sortedHits.find( northeast.rawId() );
00331   if(i != sortedHits.end() ) 
00332     rh.add8Neighbour( i->second );
00333   
00334   i = sortedHits.find( south.rawId() );
00335   if(i != sortedHits.end() ) 
00336     rh.add4Neighbour( i->second );
00337   
00338   i = sortedHits.find( southwest.rawId() );
00339   if(i != sortedHits.end() ) 
00340     rh.add8Neighbour( i->second );
00341   
00342   i = sortedHits.find( east.rawId() );
00343   if(i != sortedHits.end() ) 
00344     rh.add4Neighbour( i->second );
00345   
00346   i = sortedHits.find( southeast.rawId() );
00347   if(i != sortedHits.end() ) 
00348     rh.add8Neighbour( i->second );
00349   
00350   i = sortedHits.find( west.rawId() );
00351   if(i != sortedHits.end() ) 
00352     rh.add4Neighbour( i->second );
00353   
00354   i = sortedHits.find( northwest.rawId() );
00355   if(i != sortedHits.end() ) 
00356     rh.add8Neighbour( i->second );
00357 }
00358 
00359 
00360 // Build the array of (max)8 neighbors
00361 void 
00362 PFRecHitProducerHO::hoNeighbArray(
00363                                   const CaloSubdetectorGeometry& barrelGeom,
00364                                   const HcalTopology& barrelTopo){
00365   
00366   static const CaloDirection orderedDir[8]={SOUTHWEST,
00367                                             SOUTH,
00368                                             SOUTHEAST,
00369                                             WEST,
00370                                             EAST,
00371                                             NORTHWEST,
00372                                             NORTH,
00373                                             NORTHEAST};
00374   
00375   const unsigned nbarrel = 2160; //62000;
00376   // Barrel first. The hashed index runs from 0 to 2199 61199
00377   neighboursHO_.resize(barrelTopo.getHOSize());
00378   
00379   //std::cout << " Building the array of neighbours (barrel) " ;
00380   
00381   const std::vector<DetId>& vec(barrelGeom.getValidDetIds(DetId::Hcal,
00382                                                           HcalOuter));
00383   unsigned size=vec.size();    
00384   for(unsigned ic=0; ic<size; ++ic) 
00385     {
00386       // We get the 9 cells in a square. 
00387       std::vector<DetId> neighbours(barrelTopo.getWindow(vec[ic],3,3));
00388       unsigned nneighbours=neighbours.size();
00389       
00390       unsigned hashedindex=barrelTopo.detId2denseIdHO(vec[ic]);
00391       //           std::cout << " Cell " << ic<<" "<<vec[ic].rawId()<<" "<<HcalDetId(vec[ic]) <<" "<<hashedindex<<" "<<nneighbours<< std::endl;
00392       if(hashedindex>=nbarrel)
00393         {
00394           LogDebug("CaloGeometryTools")  << " Array overflow " << std::endl;
00395         }
00396       
00397       
00398       // If there are 9 cells, it is easy, and this order is know:
00399       //      6  7  8
00400       //      3  4  5 
00401       //      0  1  2   (0 = SOUTHWEST)
00402       
00403       if(nneighbours==9)
00404         {
00405           neighboursHO_[hashedindex].reserve(8);
00406           for(unsigned in=0;in<nneighbours;++in)
00407             {
00408               // remove the centre
00409               //    cout<<"ic "<<ic<<" "<<in<<" "<<neighbours[in].rawId()<<" "<<vec[ic].rawId()<<" "<<hashedindex<<endl; 
00410               if(neighbours[in]!=vec[ic]) 
00411                 {
00412                   neighboursHO_[hashedindex].push_back(neighbours[in]);
00413                   //            std::cout << " Neighbour " << ic<<" "<<size<<" "<<in << " " <<hashedindex<<" "<< HcalDetId(neighbours[in]) << std::endl;
00414                 }
00415             }
00416         }
00417       else
00418         {
00419           DetId central(vec[ic]);
00420           neighboursHO_[hashedindex].resize(8,DetId(0));
00421           for(unsigned idir=0;idir<8;++idir)
00422             {
00423               DetId testid=central;
00424               bool status=stdmove(testid,orderedDir[idir],
00425                                   barrelTopo, barrelGeom);
00426               if(status) neighboursHO_[hashedindex][idir]=testid;
00427             }
00428           
00429         }
00430     }
00431   
00432   //    std::cout << " done " << size <<std::endl;
00433   neighbourmapcalculated_ = true;
00434 }
00435 
00436 bool 
00437 PFRecHitProducerHO::stdsimplemove(DetId& cell, 
00438                                   const CaloDirection& dir,
00439                                   const CaloSubdetectorTopology& barrelTopo,
00440                                   const CaloSubdetectorGeometry& barrelGeom)
00441   const {
00442   
00443   std::vector<DetId> neighbours;
00444   
00445   // BARREL CASE 
00446   if(cell.subdetId()==HcalOuter) {
00447     HcalDetId hoDetId = cell;
00448     
00449     neighbours = barrelTopo.getNeighbours(hoDetId,dir);
00450     
00451     // first try to move according to the standard navigation
00452     if(neighbours.size()>0 && !neighbours[0].null()) {
00453       cell = neighbours[0];
00454       return true;
00455     }
00456     
00457     // failed.
00458     
00459     
00460   }
00461   
00462   // everything failed 
00463   cell = DetId(0);
00464   return false;
00465 }
00466 
00467 
00468 
00469 bool 
00470 PFRecHitProducerHO::stdmove(DetId& cell, 
00471                             const CaloDirection& dir,
00472                             const CaloSubdetectorTopology& barrelTopo,
00473                             const CaloSubdetectorGeometry& barrelGeom)
00474   
00475   const {
00476   
00477   
00478   bool result; 
00479   
00480   if(dir==NORTH) {
00481     result = stdsimplemove(cell,NORTH, barrelTopo, barrelGeom);
00482     return result;
00483   }
00484   else if(dir==SOUTH) {
00485     result = stdsimplemove(cell,SOUTH, barrelTopo, barrelGeom);
00486     return result;
00487   }
00488   else if(dir==EAST) {
00489     result = stdsimplemove(cell,EAST, barrelTopo, barrelGeom);
00490     return result;
00491   }
00492   else if(dir==WEST) {
00493     result = stdsimplemove(cell,WEST, barrelTopo, barrelGeom);
00494     return result;
00495   }
00496   
00497   
00498   // One has to try both paths
00499   else if(dir==NORTHEAST)
00500     {
00501       result = stdsimplemove(cell,NORTH, barrelTopo, barrelGeom);
00502       if(result)
00503         return stdsimplemove(cell,EAST, barrelTopo, barrelGeom);
00504       else
00505         {
00506           result = stdsimplemove(cell,EAST, barrelTopo, barrelGeom);
00507           if(result)
00508             return stdsimplemove(cell,NORTH, barrelTopo, barrelGeom);
00509           else
00510             return false; 
00511         }
00512     }
00513   else if(dir==NORTHWEST)
00514     {
00515       result = stdsimplemove(cell,NORTH, barrelTopo, barrelGeom );
00516       if(result)
00517         return stdsimplemove(cell,WEST, barrelTopo, barrelGeom );
00518       else
00519         {
00520           result = stdsimplemove(cell,WEST, barrelTopo, barrelGeom );
00521           if(result)
00522             return stdsimplemove(cell,NORTH, barrelTopo, barrelGeom );
00523           else
00524             return false; 
00525         }
00526     }
00527   else if(dir == SOUTHEAST)
00528     {
00529       result = stdsimplemove(cell,SOUTH, barrelTopo, barrelGeom );
00530       if(result)
00531         return stdsimplemove(cell,EAST, barrelTopo, barrelGeom );
00532       else
00533         {
00534           result = stdsimplemove(cell,EAST, barrelTopo, barrelGeom );
00535           if(result)
00536             return stdsimplemove(cell,SOUTH, barrelTopo, barrelGeom );
00537           else
00538             return false; 
00539         }
00540     }
00541   else if(dir == SOUTHWEST)
00542     {
00543       result = stdsimplemove(cell,SOUTH, barrelTopo, barrelGeom );
00544       if(result)
00545         return stdsimplemove(cell,WEST, barrelTopo, barrelGeom );
00546       else
00547         {
00548           result = stdsimplemove(cell,SOUTH, barrelTopo, barrelGeom );
00549           if(result)
00550             return stdsimplemove(cell,WEST, barrelTopo, barrelGeom );
00551           else
00552             return false; 
00553         }
00554     }
00555   cell = DetId(0);
00556   return false;
00557 }
00558 
00559 DetId PFRecHitProducerHO::move(DetId cell, 
00560                                const HcalTopology&topo,
00561                                const CaloDirection&dir ) const
00562 {  
00563   DetId originalcell = cell; 
00564   if(dir==NONE || cell==DetId(0)) return false;
00565   
00566   // Conversion CaloDirection and index in the table
00567   // CaloDirection :NONE,SOUTH,SOUTHEAST,SOUTHWEST,EAST,WEST, NORTHEAST,NORTHWEST,NORTH
00568   // Table : SOUTHWEST,SOUTH,SOUTHEAST,WEST,EAST,NORTHWEST,NORTH, NORTHEAST
00569   static const int calodirections[9]={-1,1,2,0,4,3,7,5,6};
00570   
00571   assert(neighbourmapcalculated_);
00572   
00573   DetId result = neighboursHO_[topo.detId2denseIdHO(originalcell)][calodirections[dir]];
00574   return result; 
00575 }
00576