CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/RecoMuon/L3TrackFinder/src/MuonCkfTrajectoryBuilder.cc

Go to the documentation of this file.
00001 #include "RecoMuon/L3TrackFinder/interface/MuonCkfTrajectoryBuilder.h"
00002 
00003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00004 
00005 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
00006 #include "RecoTracker/TkDetLayers/interface/GeometricSearchTracker.h"
00007 #include "TrackingTools/PatternTools/interface/TrajMeasLessEstim.h"
00008 #include "RecoTracker/MeasurementDet/interface/MeasurementTracker.h"
00009 #include "TrackingTools/MeasurementDet/interface/LayerMeasurements.h"
00010 #include "TrackingTools/KalmanUpdators/interface/Chi2MeasurementEstimatorBase.h"
00011 #include "TrackingTools/TrajectoryFiltering/interface/TrajectoryFilter.h"
00012 #include "TrackingTools/PatternTools/interface/TransverseImpactPointExtrapolator.h"
00013 #include "RecoMuon/L3TrackFinder/src/EtaPhiEstimator.h"
00014 #include <sstream>
00015 
00016 MuonCkfTrajectoryBuilder::MuonCkfTrajectoryBuilder(const edm::ParameterSet&              conf,
00017                                                    const TrajectoryStateUpdator*         updator,
00018                                                    const Propagator*                     propagatorAlong,
00019                                                    const Propagator*                     propagatorOpposite,
00020                                                    const Propagator*                     propagatorProximity,
00021                                                    const Chi2MeasurementEstimatorBase*   estimator,
00022                                                    const TransientTrackingRecHitBuilder* RecHitBuilder,
00023                                                    const MeasurementTracker*             measurementTracker,
00024                                                    const TrajectoryFilter*               filter): 
00025   CkfTrajectoryBuilder(conf,updator,propagatorAlong,propagatorOpposite,estimator,RecHitBuilder,measurementTracker,filter),
00026   theProximityPropagator(propagatorProximity)
00027 {
00028   //and something specific to me ?
00029   theUseSeedLayer = conf.getParameter<bool>("useSeedLayer");
00030   theRescaleErrorIfFail = conf.getParameter<double>("rescaleErrorIfFail");
00031   double dEta=conf.getParameter<double>("deltaEta");
00032   double dPhi=conf.getParameter<double>("deltaPhi");
00033   if (dEta>0 && dPhi>0)
00034     theEtaPhiEstimator = new EtaPhiEstimator(dEta,dPhi,theEstimator);
00035   else theEtaPhiEstimator = (Chi2MeasurementEstimatorBase*)0;
00036 
00037 
00038 }
00039 
00040 MuonCkfTrajectoryBuilder::~MuonCkfTrajectoryBuilder()
00041 {
00042   if (theEtaPhiEstimator) delete theEtaPhiEstimator;
00043 }
00044 
00045 /*
00046 std::string dumpMeasurement(const TrajectoryMeasurement & tm)
00047 {
00048   std::stringstream buffer;
00049   buffer<<"layer pointer: "<<tm.layer()<<"\n"
00050         <<"estimate: "<<tm.estimate()<<"\n"
00051         <<"forward state: \n"
00052         <<"x: "<<tm.forwardPredictedState().globalPosition()<<"\n"
00053         <<"p: "<<tm.forwardPredictedState().globalMomentum()<<"\n"
00054         <<"geomdet pointer from rechit: "<<tm.recHit()->det()<<"\n"
00055         <<"detId: "<<tm.recHit()->geographicalId().rawId();
00056   if (tm.recHit()->isValid()){
00057     buffer<<"\n hit global x: "<<tm.recHit()->globalPosition()
00058           <<"\n hit global error: "<<tm.recHit()->globalPositionError().matrix()
00059           <<"\n hit local x:"<<tm.recHit()->localPosition()
00060           <<"\n hit local error"<<tm.recHit()->localPositionError();
00061   }
00062   return buffer.str();
00063 }
00064 std::string dumpMeasurements(const std::vector<TrajectoryMeasurement> & v)
00065 {
00066   std::stringstream buffer;
00067   buffer<<v.size()<<" total measurements\n";
00068   for (std::vector<TrajectoryMeasurement>::const_iterator it = v.begin(); it!=v.end();++it){
00069     buffer<<dumpMeasurement(*it);
00070     buffer<<"\n";}
00071   return buffer.str();
00072 }
00073 */
00074 
00075 void MuonCkfTrajectoryBuilder::collectMeasurement(const DetLayer* layer,
00076                                                   const std::vector<const DetLayer*>& nl,
00077                                                   const TrajectoryStateOnSurface & currentState,
00078                                                   std::vector<TM>& result,int& invalidHits,
00079                                                   const Propagator * prop) const{
00080   for (std::vector<const DetLayer*>::const_iterator il = nl.begin();
00081        il != nl.end(); il++) {
00082 
00083     TSOS stateToUse = currentState;
00084     
00085     if (layer == (*il)){
00086       LogDebug("CkfPattern")<<" self propagating in findCompatibleMeasurements.\n from: \n"<<stateToUse;
00087       //self navigation case
00088       // go to a middle point first
00089       TransverseImpactPointExtrapolator middle;
00090       GlobalPoint center(0,0,0);
00091       stateToUse = middle.extrapolate(stateToUse, center, *prop);
00092 
00093       if (!stateToUse.isValid()) continue;
00094       LogDebug("CkfPattern")<<"to: "<<stateToUse;
00095     }
00096 
00097     std::vector<TM> tmp =
00098       theLayerMeasurements->measurements((**il),stateToUse, *prop, *theEstimator);
00099     
00100     if (tmp.size()==1 && theEtaPhiEstimator){
00101       LogDebug("CkfPattern")<<"only an invalid hit is found. trying differently";
00102       tmp = theLayerMeasurements->measurements((**il),stateToUse, *prop, *theEtaPhiEstimator);
00103     }
00104     LogDebug("CkfPattern")<<tmp.size()<<" measurements returned by LayerMeasurements";
00105     
00106     if ( !tmp.empty()) {
00107       // FIXME durty-durty-durty cleaning: never do that please !
00108       /*      for (vector<TM>::iterator it = tmp.begin(); it!=tmp.end(); ++it)
00109               {if (it->recHit()->det()==0) it=tmp.erase(it)--;}*/
00110       
00111       if ( result.empty()) result = tmp;
00112       else {
00113         // keep one dummy TM at the end, skip the others
00114         result.insert( result.end()-invalidHits, tmp.begin(), tmp.end());
00115       }
00116       invalidHits++;
00117     }
00118   }
00119   
00120   LogDebug("CkfPattern")<<"starting from:\n"
00121                         <<"x: "<<currentState.globalPosition()<<"\n"
00122                         <<"p: "<<currentState.globalMomentum()<<"\n"
00123                         <<PrintoutHelper::dumpMeasurements(result);
00124 }
00125  
00126 
00127 
00128 void 
00129 MuonCkfTrajectoryBuilder::findCompatibleMeasurements( const TempTrajectory& traj, 
00130                                                   std::vector<TrajectoryMeasurement> & result) const
00131 {
00132   int invalidHits = 0;
00133 
00134 
00135   std::vector<const DetLayer*> nl;
00136 
00137   if (traj.empty())
00138     {
00139       LogDebug("CkfPattern")<<"using JR patch for no measurement case";           
00140      //what if there are no measurement on the Trajectory
00141 
00142       //set the currentState to be the one from the trajectory seed starting point
00143       PTrajectoryStateOnDet ptod =traj.seed().startingState();
00144       DetId id(ptod.detId());
00145       const GeomDet * g = theMeasurementTracker->geomTracker()->idToDet(id);
00146       const Surface * surface=&g->surface();
00147       
00148       TrajectoryStateOnSurface currentState(trajectoryStateTransform::transientState(ptod,surface,theForwardPropagator->magneticField()));
00149 
00150       //set the next layers to be that one the state is on
00151       const DetLayer * l=theMeasurementTracker->geometricSearchTracker()->detLayer(id);
00152 
00153       if (theUseSeedLayer){
00154         {
00155           //get the measurements on the layer first
00156           LogDebug("CkfPattern")<<"using the layer of the seed first.";
00157           nl.push_back(l);
00158           collectMeasurement(l,nl,currentState,result,invalidHits,theProximityPropagator);
00159         }
00160         
00161         //if fails: try to rescale locally the state to find measurements
00162         if ((unsigned int)invalidHits==result.size() && theRescaleErrorIfFail!=1.0 && result.size()!=0)
00163           {
00164             result.clear();
00165             LogDebug("CkfPattern")<<"using a rescale by "<< theRescaleErrorIfFail <<" to find measurements.";
00166             TrajectoryStateOnSurface rescaledCurrentState = currentState;
00167             rescaledCurrentState.rescaleError(theRescaleErrorIfFail);
00168             invalidHits=0;
00169             collectMeasurement(l,nl,rescaledCurrentState,result,invalidHits,theProximityPropagator);
00170           }
00171       }
00172 
00173       //if fails: go to next layers
00174       if (result.size()==0 || (unsigned int)invalidHits==result.size())
00175         {
00176           result.clear();
00177           LogDebug("CkfPattern")<<"Need to go to next layer to get measurements";
00178           //the following will "JUMP" the first layer measurements
00179           nl = l->nextLayers(*currentState.freeState(), traj.direction());
00180           if (nl.size()==0){
00181             LogDebug("CkfPattern")<<" there was no next layer with wellInside. Use the next with no check.";
00182             //means you did not get any compatible layer on the next 1/2 tracker layer.
00183             // use the next layers with no checking
00184             nl = l->nextLayers(((traj.direction()==alongMomentum)?insideOut:outsideIn));
00185           }
00186           invalidHits=0;
00187           collectMeasurement(l,nl,currentState,result,invalidHits,theForwardPropagator);
00188         }
00189 
00190       //if fails: this is on the next layers already, try rescaling locally the state
00191       if (result.size()!=0 && (unsigned int)invalidHits==result.size() && theRescaleErrorIfFail!=1.0)
00192         {
00193           result.clear();
00194           LogDebug("CkfPattern")<<"using a rescale by "<< theRescaleErrorIfFail <<" to find measurements on next layers.";
00195           TrajectoryStateOnSurface rescaledCurrentState = currentState;
00196           rescaledCurrentState.rescaleError(theRescaleErrorIfFail);
00197           invalidHits=0;
00198           collectMeasurement(l,nl,rescaledCurrentState, result,invalidHits,theForwardPropagator);
00199         }
00200 
00201     }
00202   else //regular case
00203     {
00204 
00205       TSOS currentState( traj.lastMeasurement().updatedState());
00206 
00207       nl = traj.lastLayer()->nextLayers( *currentState.freeState(), traj.direction());
00208       if (nl.empty()){LogDebug("CkfPattern")<<" no next layers... going "<<traj.direction()<<"\n from: \n"<<currentState<<"\n from detId: "<<traj.lastMeasurement().recHit()->geographicalId().rawId(); return ;}
00209 
00210       collectMeasurement(traj.lastLayer(),nl,currentState,result,invalidHits,theForwardPropagator);
00211     }
00212 
00213 
00214   // sort the final result, keep dummy measurements at the end
00215   if ( result.size() > 1) {
00216     sort( result.begin(), result.end()-invalidHits, TrajMeasLessEstim());
00217   }
00218 
00219 #ifdef DEBUG_INVALID
00220   bool afterInvalid = false;
00221   for (std::vector<TM>::const_iterator i=result.begin();
00222        i!=result.end(); i++) {
00223     if ( ! i->recHit().isValid()) afterInvalid = true;
00224     if (afterInvalid && i->recHit().isValid()) {
00225       edm::LogError("CkfPattern") << "CkfTrajectoryBuilder error: valid hit after invalid!" ;
00226     }
00227   }
00228 #endif
00229 
00230   //analyseMeasurements( result, traj);
00231 
00232 }
00233 
00234