CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_9_patch3/src/RecoMuon/MeasurementDet/src/MuonDetLayerMeasurements.cc

Go to the documentation of this file.
00001 
00010 #include "RecoMuon/MeasurementDet/interface/MuonDetLayerMeasurements.h"
00011 
00012 #include "TrackingTools/PatternTools/interface/TrajectoryMeasurement.h" 
00013 #include "TrackingTools/DetLayers/interface/DetLayer.h"
00014 
00015 #include "TrackingTools/PatternTools/interface/TrajMeasLessEstim.h"
00016 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00017 #include "FWCore/ServiceRegistry/interface/Service.h"
00018 #include "FWCore/Services/interface/UpdaterService.h"
00019 
00020 
00021 typedef MuonTransientTrackingRecHit::MuonRecHitPointer MuonRecHitPointer;
00022 typedef MuonTransientTrackingRecHit::MuonRecHitContainer MuonRecHitContainer;
00023 
00024 
00025 
00026 MuonDetLayerMeasurements::MuonDetLayerMeasurements(edm::InputTag dtlabel, 
00027                                                    edm::InputTag csclabel, 
00028                                                    edm::InputTag rpclabel,
00029                                                    bool enableDT, bool enableCSC, bool enableRPC): 
00030   theDTRecHitLabel(dtlabel),
00031   theCSCRecHitLabel(csclabel),
00032   theRPCRecHitLabel(rpclabel),
00033   enableDTMeasurement(enableDT),
00034   enableCSCMeasurement(enableCSC),
00035   enableRPCMeasurement(enableRPC),
00036   theDTRecHits(),
00037   theCSCRecHits(),
00038   theRPCRecHits(),
00039   theDTEventID(),
00040   theCSCEventID(),
00041   theRPCEventID(),
00042   theEvent(0){
00043   static int procInstance(0);
00044   std::ostringstream sDT;
00045   sDT<<"MuonDetLayerMeasurements::checkDTRecHits::" << procInstance;
00046   theDTCheckName = sDT.str();
00047   std::ostringstream sRPC;
00048   sRPC<<"MuonDetLayerMeasurements::checkRPCRecHits::" << procInstance;
00049   theRPCCheckName = sRPC.str();
00050   std::ostringstream sCSC;
00051   sCSC<<"MuonDetLayerMeasurements::checkCSCRecHits::" << procInstance;
00052   theCSCCheckName = sCSC.str();
00053   procInstance++;
00054 }
00055 
00056 MuonDetLayerMeasurements::~MuonDetLayerMeasurements(){}
00057 
00058 MuonRecHitContainer MuonDetLayerMeasurements::recHits(const GeomDet* geomDet, 
00059                                                       const edm::Event& iEvent)
00060 {
00061   DetId geoId = geomDet->geographicalId();
00062   theEvent = &iEvent;
00063   MuonRecHitContainer result;
00064 
00065   if (geoId.subdetId()  == MuonSubdetId::DT) {
00066     if(enableDTMeasurement) 
00067     {
00068       checkDTRecHits();
00069     
00070       // Create the ChamberId
00071       DTChamberId chamberId(geoId.rawId());
00072       // LogTrace("Muon|RecoMuon|MuonDetLayerMeasurements") << "(DT): "<<chamberId<<std::endl;
00073     
00074       // Get the DT-Segment which relies on this chamber
00075       DTRecSegment4DCollection::range range = theDTRecHits->get(chamberId);
00076     
00077       // Create the MuonTransientTrackingRechit
00078       for (DTRecSegment4DCollection::const_iterator rechit = range.first; 
00079            rechit!=range.second;++rechit)
00080         result.push_back(MuonTransientTrackingRecHit::specificBuild(geomDet,&*rechit));
00081     }
00082   }
00083   
00084   else if (geoId.subdetId()  == MuonSubdetId::CSC) {
00085     if(enableCSCMeasurement)
00086     {
00087       checkCSCRecHits();
00088 
00089       // Create the chamber Id
00090       CSCDetId chamberId(geoId.rawId());
00091       //    LogTrace("Muon|RecoMuon|MuonDetLayerMeasurements") << "(CSC): "<<chamberId<<std::endl;
00092 
00093       // Get the CSC-Segment which relies on this chamber
00094       CSCSegmentCollection::range range = theCSCRecHits->get(chamberId);
00095     
00096       // Create the MuonTransientTrackingRecHit
00097       for (CSCSegmentCollection::const_iterator rechit = range.first; 
00098            rechit!=range.second; ++rechit)
00099         result.push_back(MuonTransientTrackingRecHit::specificBuild(geomDet,&*rechit)); 
00100     }
00101   }
00102   
00103   else if (geoId.subdetId()  == MuonSubdetId::RPC) {
00104     if(enableRPCMeasurement)
00105     {
00106       checkRPCRecHits(); 
00107 
00108       // Create the chamber Id
00109       RPCDetId chamberId(geoId.rawId());
00110       // LogTrace("Muon|RecoMuon|MuonDetLayerMeasurements") << "(RPC): "<<chamberId<<std::endl;
00111     
00112       // Get the RPC-Segment which relies on this chamber
00113       RPCRecHitCollection::range range = theRPCRecHits->get(chamberId);
00114     
00115       // Create the MuonTransientTrackingRecHit
00116       for (RPCRecHitCollection::const_iterator rechit = range.first; 
00117            rechit!=range.second; ++rechit)
00118         result.push_back(MuonTransientTrackingRecHit::specificBuild(geomDet,&*rechit));
00119     }
00120   }
00121   else {
00122     // wrong type
00123     throw cms::Exception("MuonDetLayerMeasurements") << "The DetLayer with det " << geoId.det() << " subdet " << geoId.subdetId() << " is not a valid Muon DetLayer. ";
00124   }
00125   return result;
00126 }
00127 
00128 
00129 void MuonDetLayerMeasurements::checkDTRecHits()
00130 {
00131   checkEvent();
00132   if (!edm::Service<UpdaterService>()->checkOnce(theDTCheckName)) return;
00133 
00134   {
00135     theDTEventID = theEvent->id();
00136     theEvent->getByLabel(theDTRecHitLabel, theDTRecHits);
00137   }
00138   if(!theDTRecHits.isValid())
00139   {
00140     throw cms::Exception("MuonDetLayerMeasurements") << "Cannot get DT RecHits";
00141   }
00142 }
00143 
00144 
00145 void MuonDetLayerMeasurements::checkCSCRecHits()
00146 {
00147   checkEvent();
00148   if (!edm::Service<UpdaterService>()->checkOnce(theCSCCheckName)) return;
00149 
00150   {
00151     theCSCEventID = theEvent->id();
00152     theEvent->getByLabel(theCSCRecHitLabel, theCSCRecHits);
00153   }
00154   if(!theCSCRecHits.isValid())
00155   {
00156     throw cms::Exception("MuonDetLayerMeasurements") << "Cannot get CSC RecHits";
00157   }
00158 }
00159 
00160 
00161 void MuonDetLayerMeasurements::checkRPCRecHits()
00162 {
00163   checkEvent();
00164   if (!edm::Service<UpdaterService>()->checkOnce(theRPCCheckName)) return;
00165 
00166   {
00167     theRPCEventID = theEvent->id();
00168     theEvent->getByLabel(theRPCRecHitLabel, theRPCRecHits);
00169   }
00170   if(!theRPCRecHits.isValid())
00171   {
00172     throw cms::Exception("MuonDetLayerMeasurements") << "Cannot get RPC RecHits";
00173   }
00174 }
00175 
00176 
00178 MeasurementContainer
00179 MuonDetLayerMeasurements::measurements( const DetLayer* layer,
00180                                         const TrajectoryStateOnSurface& startingState,
00181                                         const Propagator& prop,
00182                                         const MeasurementEstimator& est) {
00183   checkEvent();
00184   return measurements(layer, startingState, prop, est, *theEvent);
00185 }
00186 
00187 
00188 MeasurementContainer
00189 MuonDetLayerMeasurements::measurements(const DetLayer* layer,
00190                                        const TrajectoryStateOnSurface& startingState,
00191                                        const Propagator& prop,
00192                                        const MeasurementEstimator& est,
00193                                        const edm::Event& iEvent) {
00194   
00195   MeasurementContainer result;
00196   
00197   std::vector<DetWithState> dss = layer->compatibleDets(startingState, prop, est);
00198   LogTrace("RecoMuon")<<"compatibleDets: "<<dss.size()<<std::endl;
00199   
00200   for(std::vector<DetWithState>::const_iterator detWithStateItr = dss.begin();
00201       detWithStateItr != dss.end(); ++detWithStateItr){
00202 
00203     MeasurementContainer detMeasurements 
00204       = measurements(layer, detWithStateItr->first, 
00205                      detWithStateItr->second, est, iEvent);
00206     result.insert(result.end(), detMeasurements.begin(), detMeasurements.end());
00207   }
00208   
00209   if (!result.empty()) sort( result.begin(), result.end(), TrajMeasLessEstim());
00210   
00211   return result;
00212 }
00213 
00214 
00215 MeasurementContainer
00216 MuonDetLayerMeasurements::measurements( const DetLayer* layer,
00217                                         const GeomDet* det,
00218                                         const TrajectoryStateOnSurface& stateOnDet,
00219                                         const MeasurementEstimator& est,
00220                                         const edm::Event& iEvent) {
00221   MeasurementContainer result;
00222   
00223   // Get the Segments which relies on the GeomDet given by compatibleDets
00224   MuonRecHitContainer muonRecHits = recHits(det, iEvent);
00225   
00226   // Create the Trajectory Measurement
00227   for(MuonRecHitContainer::const_iterator rechit = muonRecHits.begin();
00228       rechit != muonRecHits.end(); ++rechit) {
00229 
00230     MeasurementEstimator::HitReturnType estimate = est.estimate(stateOnDet,**rechit);
00231     LogTrace("RecoMuon")<<"Dimension: "<<(*rechit)->dimension()
00232                         <<" Chi2: "<<estimate.second<<std::endl;
00233     if (estimate.first) {
00234       result.push_back(TrajectoryMeasurement(stateOnDet, rechit->get(),
00235                                              estimate.second,layer));
00236     }
00237   }
00238 
00239   if (!result.empty()) sort( result.begin(), result.end(), TrajMeasLessEstim());
00240    
00241   return result;
00242 }
00243 
00244 
00245 
00246 MeasurementContainer
00247 MuonDetLayerMeasurements::fastMeasurements( const DetLayer* layer,
00248                                             const TrajectoryStateOnSurface& theStateOnDet,
00249                                             const TrajectoryStateOnSurface& startingState,
00250                                             const Propagator& prop,
00251                                             const MeasurementEstimator& est,
00252                                             const edm::Event& iEvent) {
00253   MeasurementContainer result;
00254   MuonRecHitContainer rhs = recHits(layer, iEvent);
00255   for (MuonRecHitContainer::const_iterator irh = rhs.begin(); irh!=rhs.end(); irh++) {
00256     MeasurementEstimator::HitReturnType estimate = est.estimate(theStateOnDet, (**irh));
00257     if (estimate.first)
00258     {
00259       result.push_back(TrajectoryMeasurement(theStateOnDet,(*irh).get(),
00260                                              estimate.second,layer));
00261     }
00262   }
00263 
00264   if (!result.empty()) {
00265     sort( result.begin(), result.end(), TrajMeasLessEstim());
00266   }
00267 
00268   return result;
00269 }
00270 
00272 MeasurementContainer
00273 MuonDetLayerMeasurements::fastMeasurements(const DetLayer* layer,
00274                                            const TrajectoryStateOnSurface& theStateOnDet,
00275                                            const TrajectoryStateOnSurface& startingState,
00276                                            const Propagator& prop,
00277                                            const MeasurementEstimator& est) {
00278   checkEvent();
00279   return fastMeasurements(layer, theStateOnDet, startingState, prop, est, *theEvent); 
00280 }
00281 
00282 
00283 std::vector<TrajectoryMeasurementGroup>
00284 MuonDetLayerMeasurements::groupedMeasurements(const DetLayer* layer,
00285                                               const TrajectoryStateOnSurface& startingState,
00286                                               const Propagator& prop,
00287                                               const MeasurementEstimator& est) {
00288   checkEvent();
00289   return groupedMeasurements(layer, startingState, prop,  est, *theEvent);
00290 }
00291 
00292 
00293 std::vector<TrajectoryMeasurementGroup>
00294 MuonDetLayerMeasurements::groupedMeasurements(const DetLayer* layer,
00295                                               const TrajectoryStateOnSurface& startingState,
00296                                               const Propagator& prop,
00297                                               const MeasurementEstimator& est,
00298                                               const edm::Event& iEvent) {
00299   
00300   std::vector<TrajectoryMeasurementGroup> result;
00301   // if we want to use the concept of InvalidRecHits,
00302   // we can reuse LayerMeasurements from TrackingTools/MeasurementDet
00303   std::vector<DetGroup> groups(layer->groupedCompatibleDets(startingState, prop, est));
00304 
00305   // this should be fixed either in RecoMuon/MeasurementDet/MuonDetLayerMeasurements or
00306   // RecoMuon/DetLayers/MuRingForwardDoubleLayer
00307   // and removed the reverse operation in StandAloneMuonFilter::findBestMeasurements
00308 
00309   for (std::vector<DetGroup>::const_iterator grp=groups.begin(); grp!=groups.end(); ++grp) {
00310     
00311     std::vector<TrajectoryMeasurement> groupMeasurements;
00312     for (DetGroup::const_iterator detAndStateItr=grp->begin();
00313          detAndStateItr !=grp->end(); ++detAndStateItr) {
00314 
00315       std::vector<TrajectoryMeasurement> detMeasurements 
00316         = measurements(layer, detAndStateItr->det(), detAndStateItr->trajectoryState(), est, iEvent);
00317       groupMeasurements.insert(groupMeasurements.end(), detMeasurements.begin(), detMeasurements.end());
00318     }
00319     
00320     if (!groupMeasurements.empty()) 
00321       std::sort( groupMeasurements.begin(), groupMeasurements.end(), TrajMeasLessEstim());  
00322     
00323     result.push_back(TrajectoryMeasurementGroup(groupMeasurements, *grp));
00324   }
00325 
00326   return result;
00327 }
00328 
00330 void MuonDetLayerMeasurements::setEvent(const edm::Event& event) {
00331   theEvent = &event;
00332 }
00333 
00334 
00335 void MuonDetLayerMeasurements::checkEvent() const {
00336   if(!theEvent)
00337     throw cms::Exception("MuonDetLayerMeasurements") << "The event has not been set";
00338 }
00339 
00340 MuonRecHitContainer MuonDetLayerMeasurements::recHits(const DetLayer* layer, 
00341                                                       const edm::Event& iEvent) {
00342   MuonRecHitContainer rhs;
00343   
00344   std::vector <const GeomDet*> gds = layer->basicComponents();
00345 
00346   for (std::vector<const GeomDet*>::const_iterator igd = gds.begin(); 
00347        igd != gds.end(); igd++) {
00348     MuonRecHitContainer detHits = recHits(*igd, iEvent);
00349     rhs.insert(rhs.end(), detHits.begin(), detHits.end());
00350   }
00351   return rhs;
00352 }
00353 
00354 MuonRecHitContainer MuonDetLayerMeasurements::recHits(const DetLayer* layer) 
00355 {
00356   checkEvent();
00357   return recHits(layer, *theEvent);
00358 }
00359