CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10/src/FastSimulation/Muons/plugins/FastTSGFromPropagation.cc

Go to the documentation of this file.
00001 #include "FastSimulation/Muons/plugins/FastTSGFromPropagation.h"
00002 
00010 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
00011 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
00012 #include "TrackingTools/PatternTools/interface/Trajectory.h"
00013 #include "TrackingTools/MeasurementDet/interface/LayerMeasurements.h"
00014 #include "TrackingTools/MeasurementDet/interface/MeasurementDet.h"
00015 #include "TrackingTools/KalmanUpdators/interface/Chi2MeasurementEstimator.h"
00016 #include "TrackingTools/GeomPropagators/interface/Propagator.h"
00017 #include "TrackingTools/GeomPropagators/interface/StateOnTrackerBound.h"
00018 
00019 #include "RecoTracker/Record/interface/TrackerRecoGeometryRecord.h"
00020 #include "RecoTracker/Record/interface/CkfComponentsRecord.h"
00021 // #include "RecoTracker/MeasurementDet/interface/TkStripMeasurementDet.h"
00022 #include "RecoTracker/MeasurementDet/interface/MeasurementTracker.h"
00023 #include "RecoTracker/TkDetLayers/interface/GeometricSearchTracker.h"
00024 
00025 #include "RecoMuon/GlobalTrackingTools/interface/DirectTrackerNavigation.h"
00026 #include "TrackingTools/KalmanUpdators/interface/KFUpdator.h"
00027 
00028 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00029 #include "FWCore/ServiceRegistry/interface/Service.h"
00030 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00031 
00032 #include "SimDataFormats/Track/interface/SimTrackContainer.h"
00033 #include "DataFormats/TrackerRecHit2D/interface/SiTrackerGSRecHit2DCollection.h"
00034 #include "DataFormats/TrackerRecHit2D/interface/SiTrackerGSMatchedRecHit2DCollection.h"
00035 #include "FastSimulation/Tracking/interface/TrackerRecHit.h"
00036 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00037 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00038 #include "Geometry/CommonDetUnit/interface/GeomDetUnit.h"
00039 #include "FastSimulation/Tracking/plugins/TrajectorySeedProducer.h"
00040 #include "TrackingTools/TransientTrackingRecHit/interface/GenericTransientTrackingRecHit.h"
00041 #include "TrackingTools/Records/interface/TransientRecHitRecord.h"
00042 
00043 
00044 using namespace std;
00045 
00046 
00047 FastTSGFromPropagation::FastTSGFromPropagation(const edm::ParameterSet & iConfig) :theTkLayerMeasurements (0), theTracker(0), theNavigation(0), theService(0), theEstimator(0),  theSigmaZ(0), theConfig (iConfig)
00048 {
00049   theCategory = "FastSimulation|Muons||FastTSGFromPropagation";
00050 
00051 }
00052 
00053 FastTSGFromPropagation::FastTSGFromPropagation(const edm::ParameterSet & iConfig, const MuonServiceProxy* service) : theTkLayerMeasurements (0), theTracker(0), theNavigation(0), theService(service),theUpdator(0), theEstimator(0), theSigmaZ(0), theConfig (iConfig)
00054 {
00055   theCategory = "FastSimulation|Muons|FastTSGFromPropagation";
00056 }
00057 
00058 FastTSGFromPropagation::~FastTSGFromPropagation()
00059 {
00060 
00061   LogTrace(theCategory) << " FastTSGFromPropagation dtor called ";
00062   if ( theNavigation ) delete theNavigation;
00063   if ( theUpdator ) delete theUpdator;
00064   if ( theEstimator ) delete theEstimator;
00065   if ( theTkLayerMeasurements ) delete theTkLayerMeasurements;
00066   if ( theErrorMatrixAdjuster ) delete theErrorMatrixAdjuster;
00067 
00068 }
00069 
00070 void FastTSGFromPropagation::trackerSeeds(const TrackCand& staMuon, const TrackingRegion& region, std::vector<TrajectorySeed> & result) {
00071 
00072   if ( theResetMethod == "discrete" ) getRescalingFactor(staMuon);
00073 
00074   TrajectoryStateOnSurface staState = outerTkState(staMuon);
00075 
00076   if ( !staState.isValid() ) { 
00077     LogTrace(theCategory) << "Error: initial state from L2 muon is invalid.";
00078     return;
00079   }
00080 
00081   LogTrace(theCategory) << "begin of trackerSeed:\n staState pos: "<<staState.globalPosition()
00082                         << " mom: "<<staState.globalMomentum() 
00083                         <<"pos eta: "<<staState.globalPosition().eta()
00084                         <<"mom eta: "<<staState.globalMomentum().eta();
00085 
00086 
00087   std::vector<const DetLayer*> nls = theNavigation->compatibleLayers(*(staState.freeState()), oppositeToMomentum);
00088 
00089   LogTrace(theCategory) << " compatible layers: "<<nls.size();
00090 
00091   if ( nls.empty() ) return;
00092 
00093   int ndesLayer = 0;
00094 
00095   bool usePredictedState = false;
00096 
00097   if ( theUpdateStateFlag ) { //use updated states
00098      std::vector<TrajectoryMeasurement> alltm; 
00099 
00100      for (std::vector<const DetLayer*>::const_iterator inl = nls.begin();
00101          inl != nls.end(); inl++, ndesLayer++ ) {
00102          if ( (*inl == 0) ) break;
00103 //         if ( (inl != nls.end()-1 ) && ( (*inl)->subDetector() == GeomDetEnumerators::TEC ) && ( (*(inl+1))->subDetector() == GeomDetEnumerators::TOB ) ) continue; 
00104          alltm = findMeasurements_new(*inl, staState);
00105          if ( (!alltm.empty()) ) {
00106             LogTrace(theCategory) << "final compatible layer: "<<ndesLayer;
00107             break;
00108          }
00109      }
00110 
00111      if ( alltm.empty() ) {
00112         LogTrace(theCategory) << " NO Measurements Found: eta: "<<staState.globalPosition().eta() <<"pt "<<staState.globalMomentum().perp();
00113         usePredictedState = true;
00114      } else {
00115        LogTrace(theCategory) << " Measurements for seeds: "<<alltm.size();
00116        std::stable_sort(alltm.begin(),alltm.end(),increasingEstimate()); 
00117        if ( alltm.size() > 5 ) alltm.erase(alltm.begin() + 5, alltm.end());
00118 
00119        const edm::SimTrackContainer* simTracks = &(*theSimTracks);
00120        const std::vector<unsigned> theSimTrackIds = theGSRecHits->ids();
00121        TrackerRecHit theSeedHits;
00122        std::vector<TrackerRecHit> outerHits;
00123 
00124        //std::vector<TrajectorySeed>  tmpTS;
00125        bool isMatch = false;
00126        for (std::vector<TrajectoryMeasurement>::const_iterator itm = alltm.begin(); itm != alltm.end(); itm++) {
00127            const TrajectoryStateOnSurface seedState = itm->predictedState();
00128            double preY = seedState.globalPosition().y();
00129 
00130            // Check SimTrack
00131            TrackingRecHit* aTrackingRecHit;
00132            FreeTrajectoryState simtrack_trackerstate;
00133            for( unsigned tkId=0;  tkId != theSimTrackIds.size(); ++tkId ) {
00134                const SimTrack & simtrack = (*simTracks)[theSimTrackIds[tkId]];
00135                SiTrackerGSMatchedRecHit2DCollection::range theRecHitRange = theGSRecHits->get(theSimTrackIds[tkId]);
00136                SiTrackerGSMatchedRecHit2DCollection::const_iterator theRecHitRangeIteratorBegin = theRecHitRange.first;
00137                SiTrackerGSMatchedRecHit2DCollection::const_iterator theRecHitRangeIteratorEnd   = theRecHitRange.second;
00138                SiTrackerGSMatchedRecHit2DCollection::const_iterator iterRecHit;
00139 
00140                GlobalPoint position(simtrack.trackerSurfacePosition().x(),
00141                                     simtrack.trackerSurfacePosition().y(),
00142                                     simtrack.trackerSurfacePosition().z());
00143                GlobalVector momentum(simtrack.trackerSurfaceMomentum().x(),
00144                                      simtrack.trackerSurfaceMomentum().y(),
00145                                      simtrack.trackerSurfaceMomentum().z());
00146                int charge = (int)simtrack.charge();
00147                GlobalTrajectoryParameters glb_parameters(position, momentum, charge, &*theService->magneticField().product());
00148                simtrack_trackerstate = FreeTrajectoryState(glb_parameters);
00149 
00150                unsigned int outerId = 0;
00151                for( iterRecHit = theRecHitRangeIteratorBegin; iterRecHit != theRecHitRangeIteratorEnd; ++iterRecHit) {
00152                    theSeedHits = TrackerRecHit(&(*iterRecHit), theGeometry);
00153                    unsigned int id = theSeedHits.hit()->geographicalId().rawId();
00154                    if( preY < 0 ) {
00155                        if( id > outerId ) outerId = id;
00156                    }
00157                    else {
00158                        if( id > outerId ) outerId = id;
00159                    }
00160                }
00161                for( iterRecHit = theRecHitRangeIteratorBegin; iterRecHit != theRecHitRangeIteratorEnd; ++iterRecHit) {
00162                    theSeedHits = TrackerRecHit(&(*iterRecHit), theGeometry);
00163                    if( itm->recHit()->hit()->geographicalId().rawId() == theSeedHits.hit()->geographicalId().rawId() ) {
00164                        aTrackingRecHit = theSeedHits.hit()->clone();
00165                        TransientTrackingRecHit::ConstRecHitPointer recHit = theTTRHBuilder->build(aTrackingRecHit);
00166                        if( !recHit ) continue;
00167                        TrajectoryStateOnSurface updatedTSOS = updator()->update(seedState, *(recHit));
00168                        if( updatedTSOS.isValid() && passSelection(updatedTSOS) ) {
00169                            edm::OwnVector<TrackingRecHit> container;
00170                            container.push_back(recHit->hit()->clone());
00171                            TrajectorySeed ts = createSeed(updatedTSOS, container, recHit->geographicalId());
00172                            // check direction
00173                            const BasicTrajectorySeed* aSeed = &ts;
00174                            PTrajectoryStateOnDet PTSOD = aSeed->startingState();
00175                            
00176                            const GeomDet *g = theGeometry->idToDet(PTSOD.detId());
00177                            TrajectoryStateOnSurface tsos = trajectoryStateTransform::transientState(PTSOD, &(g->surface()),  &*theService->magneticField().product());
00178                            if( tsos.globalMomentum().basicVector()*seedState.globalMomentum().basicVector() < 0. ) continue;
00179                            result.push_back(ts);
00180                            isMatch = true;
00181                        }
00182                    }
00183                }
00184            }
00185        }
00186        if( !isMatch ) {
00187          // if there is no hits w.r.t. TM, find outermost hit
00188          for (std::vector<TrajectoryMeasurement>::const_iterator itm = alltm.begin(); itm != alltm.end(); itm++) {
00189            const TrajectoryStateOnSurface seedState = itm->predictedState();
00190            double preY = seedState.globalPosition().y();
00191 
00192            // Check SimTrack
00193            TrackingRecHit* aTrackingRecHit;
00194            FreeTrajectoryState simtrack_trackerstate;
00195            for( unsigned tkId=0;  tkId != theSimTrackIds.size(); ++tkId ) {
00196                const SimTrack & simtrack = (*simTracks)[theSimTrackIds[tkId]];
00197                SiTrackerGSMatchedRecHit2DCollection::range theRecHitRange = theGSRecHits->get(theSimTrackIds[tkId]);
00198                SiTrackerGSMatchedRecHit2DCollection::const_iterator theRecHitRangeIteratorBegin = theRecHitRange.first;
00199                SiTrackerGSMatchedRecHit2DCollection::const_iterator theRecHitRangeIteratorEnd   = theRecHitRange.second;
00200                SiTrackerGSMatchedRecHit2DCollection::const_iterator iterRecHit;
00201 
00202                GlobalPoint position(simtrack.trackerSurfacePosition().x(),
00203                                     simtrack.trackerSurfacePosition().y(),
00204                                     simtrack.trackerSurfacePosition().z());
00205                GlobalVector momentum(simtrack.trackerSurfaceMomentum().x(),
00206                                      simtrack.trackerSurfaceMomentum().y(),
00207                                      simtrack.trackerSurfaceMomentum().z());
00208                int charge = (int)simtrack.charge();
00209                GlobalTrajectoryParameters glb_parameters(position, momentum, charge, &*theService->magneticField().product());
00210                simtrack_trackerstate = FreeTrajectoryState(glb_parameters);
00211 
00212                unsigned int outerId = 0;
00213                for( iterRecHit = theRecHitRangeIteratorBegin; iterRecHit != theRecHitRangeIteratorEnd; ++iterRecHit) {
00214                    theSeedHits = TrackerRecHit(&(*iterRecHit), theGeometry);
00215                    unsigned int id = theSeedHits.hit()->geographicalId().rawId();
00216                    if( preY < 0 ) {
00217                        if( id > outerId ) outerId = id;
00218                    }
00219                    else {
00220                        if( id > outerId ) outerId = id;
00221                    }
00222                }
00223                for( iterRecHit = theRecHitRangeIteratorBegin; iterRecHit != theRecHitRangeIteratorEnd; ++iterRecHit) {
00224                    theSeedHits = TrackerRecHit(&(*iterRecHit), theGeometry);
00225                    if( outerId == theSeedHits.hit()->geographicalId().rawId() ) {
00226                        aTrackingRecHit = theSeedHits.hit()->clone();
00227                        TransientTrackingRecHit::ConstRecHitPointer recHit = theTTRHBuilder->build(aTrackingRecHit);
00228                        if( !recHit ) continue;
00229                        TrajectoryStateOnSurface updatedTSOS = updator()->update(seedState, *(recHit));
00230                        if( updatedTSOS.isValid() && passSelection(updatedTSOS) ) {
00231                            edm::OwnVector<TrackingRecHit> container;
00232                            container.push_back(recHit->hit()->clone());
00233                            TrajectorySeed ts = createSeed(updatedTSOS, container, recHit->geographicalId());
00234                            // check direction
00235                            const BasicTrajectorySeed* aSeed = &ts;
00236                            PTrajectoryStateOnDet PTSOD = aSeed->startingState();
00237                            
00238                            const GeomDet *g = theGeometry->idToDet(PTSOD.detId());
00239                            TrajectoryStateOnSurface tsos = trajectoryStateTransform::transientState(PTSOD, &(g->surface()),  &*theService->magneticField().product());
00240                            if( tsos.globalMomentum().basicVector()*seedState.globalMomentum().basicVector() < 0. ) continue;
00241                            result.push_back(ts);
00242                        }
00243                    }
00244                }
00245            }
00246          }
00247        }
00248 
00249        /*
00250        for( unsigned ir = 0; ir < tmpTS.size(); ir++ ) {
00251          const BasicTrajectorySeed* aSeed = &((tmpTS)[ir]);
00252          PTrajectoryStateOnDet PTSOD = aSeed->startingState();
00253          
00254          DetId seedDetId(PTSOD.detId());
00255          const GeomDet * g = theGeometry->idToDet(seedDetId);
00256          TrajectoryStateOnSurface tsos = trajectoryStateTransform::transientState(PTSOD, &(g->surface()),  &*theService->magneticField().product());
00257                  cout << "tsos3 = " << tsos.globalMomentum() << endl;
00258          if( _index == ir ) {
00259                  cout << "tsos4 = " << tsos.globalMomentum() << endl;
00260              result.push_back(tmpTS[ir]);
00261          }
00262        }
00263        */
00264        LogTrace(theCategory) << "result: "<<result.size();
00265        return;
00266      }
00267   }
00268   
00269   if ( !theUpdateStateFlag || usePredictedState ) { //use predicted states
00270      LogTrace(theCategory) << "use predicted state: ";
00271      for (std::vector<const DetLayer*>::const_iterator inl = nls.begin();
00272          inl != nls.end(); inl++ ) {
00273 
00274          if ( !result.empty() || *inl == 0 ) {
00275             break;
00276          }
00277          std::vector<DetLayer::DetWithState> compatDets = (*inl)->compatibleDets(staState, *propagator(), *estimator());
00278          LogTrace(theCategory) << " compatDets "<<compatDets.size();
00279          if ( compatDets.empty() ) continue;
00280          TrajectorySeed ts = createSeed(compatDets.front().second, compatDets.front().first->geographicalId());
00281          result.push_back(ts);
00282 
00283      }
00284      LogTrace(theCategory) << "result: "<<result.size();
00285      return;
00286   } 
00287   return;
00288 }
00289 
00290 void FastTSGFromPropagation::init(const MuonServiceProxy* service) {
00291 
00292   theMaxChi2 = theConfig.getParameter<double>("MaxChi2");
00293 
00294   theFixedErrorRescaling = theConfig.getParameter<double>("ErrorRescaling");
00295 
00296   theFlexErrorRescaling = 1.0;
00297 
00298   theResetMethod = theConfig.getParameter<std::string>("ResetMethod");
00299 
00300   if (theResetMethod != "discrete" && theResetMethod != "fixed" && theResetMethod != "matrix"  ) {
00301     edm::LogError("FastTSGFromPropagation") 
00302       <<"Wrong error rescaling method: "<<theResetMethod <<"\n"
00303       <<"Possible choices are: discrete, fixed, matrix.\n"
00304       <<"Use discrete method" <<std::endl;
00305     theResetMethod = "discrete"; 
00306   }
00307 
00308   theEstimator = new Chi2MeasurementEstimator(theMaxChi2);
00309 
00310   theCacheId_MT = 0;
00311 
00312   theCacheId_TG = 0;
00313 
00314   thePropagatorName = theConfig.getParameter<std::string>("Propagator");
00315 
00316   theService = service;
00317 
00318   theUseVertexStateFlag = theConfig.getParameter<bool>("UseVertexState");
00319 
00320   theUpdateStateFlag = theConfig.getParameter<bool>("UpdateState");
00321 
00322   theSelectStateFlag = theConfig.getParameter<bool>("SelectState");
00323 
00324   theSimTrackCollectionLabel = theConfig.getParameter<edm::InputTag>("SimTrackCollectionLabel");
00325   theHitProducer = theConfig.getParameter<edm::InputTag>("HitProducer");
00326 
00327   theUpdator = new KFUpdator();
00328 
00329   theSigmaZ = theConfig.getParameter<double>("SigmaZ");
00330 
00331   edm::ParameterSet errorMatrixPset = theConfig.getParameter<edm::ParameterSet>("errorMatrixPset");
00332   if ( theResetMethod == "matrix" && !errorMatrixPset.empty()){
00333     theAdjustAtIp = errorMatrixPset.getParameter<bool>("atIP");
00334     theErrorMatrixAdjuster = new MuonErrorMatrix(errorMatrixPset);
00335   } else {
00336     theAdjustAtIp =false;
00337     theErrorMatrixAdjuster=0;
00338   }
00339 
00340   theService->eventSetup().get<TrackerRecoGeometryRecord>().get(theTracker); 
00341   theNavigation = new DirectTrackerNavigation(theTracker);
00342 
00343   edm::ESHandle<TrackerGeometry>        geometry;
00344   theService->eventSetup().get<TrackerDigiGeometryRecord>().get(geometry);
00345   theGeometry = &(*geometry);
00346 
00347   theService->eventSetup().get<TransientRecHitRecord>().get("WithTrackAngle", theTTRHBuilder);
00348 
00349 }
00350 
00351 void FastTSGFromPropagation::setEvent(const edm::Event& iEvent) {
00352 
00353   bool measTrackerChanged = false;
00354 
00355   iEvent.getByType(theBeamSpot);
00356   
00357   // retrieve the MC truth (SimTracks)
00358   iEvent.getByLabel(theSimTrackCollectionLabel, theSimTracks);
00359   iEvent.getByLabel(theHitProducer, theGSRecHits);
00360 
00361 
00362   unsigned long long newCacheId_MT = theService->eventSetup().get<CkfComponentsRecord>().cacheIdentifier();
00363 
00364   if ( theUpdateStateFlag && newCacheId_MT != theCacheId_MT ) {
00365     LogTrace(theCategory) << "Measurment Tracker Geometry changed!";
00366     theCacheId_MT = newCacheId_MT;
00367     theService->eventSetup().get<CkfComponentsRecord>().get(theMeasTracker);
00368     measTrackerChanged = true;
00369   }
00370 
00371   //if ( theUpdateStateFlag ) theMeasTracker->update(iEvent);
00372 
00373   if ( measTrackerChanged && (&*theMeasTracker) ) {
00374      if ( theTkLayerMeasurements ) delete theTkLayerMeasurements;
00375      theTkLayerMeasurements = new LayerMeasurements(&*theMeasTracker);
00376   }
00377 
00378   bool trackerGeomChanged = false;
00379 
00380   unsigned long long newCacheId_TG = theService->eventSetup().get<TrackerRecoGeometryRecord>().cacheIdentifier();
00381 
00382   if ( newCacheId_TG != theCacheId_TG ) {
00383     LogTrace(theCategory) << "Tracker Reco Geometry changed!";
00384     theCacheId_TG = newCacheId_TG;
00385     theService->eventSetup().get<TrackerRecoGeometryRecord>().get(theTracker);
00386     trackerGeomChanged = true;
00387   }
00388 
00389   if ( trackerGeomChanged && (&*theTracker) ) {
00390     if ( theNavigation ) delete theNavigation;
00391     theNavigation = new DirectTrackerNavigation(theTracker);
00392   }
00393 }
00394 
00395 TrajectoryStateOnSurface FastTSGFromPropagation::innerState(const TrackCand& staMuon) const {
00396 
00397   TrajectoryStateOnSurface innerTS;
00398 
00399   if ( staMuon.first && staMuon.first->isValid() ) {
00400     if (staMuon.first->direction() == alongMomentum) {
00401       innerTS = staMuon.first->firstMeasurement().updatedState();
00402     } 
00403     else if (staMuon.first->direction() == oppositeToMomentum) { 
00404       innerTS = staMuon.first->lastMeasurement().updatedState();
00405     }
00406   } else {
00407     innerTS = trajectoryStateTransform::innerStateOnSurface(*(staMuon.second),*theService->trackingGeometry(), &*theService->magneticField());
00408   }
00409   //rescale the error
00410   adjust(innerTS);
00411 
00412   return  innerTS;
00413 
00414 }
00415 
00416 TrajectoryStateOnSurface FastTSGFromPropagation::outerTkState(const TrackCand& staMuon) const {
00417 
00418   TrajectoryStateOnSurface result;
00419 
00420   if ( theUseVertexStateFlag && staMuon.second->pt() > 1.0 ) {
00421     FreeTrajectoryState iniState = trajectoryStateTransform::initialFreeState(*(staMuon.second), &*theService->magneticField());
00422     //rescale the error at IP
00423     adjust(iniState); 
00424 
00425     StateOnTrackerBound fromInside(&*(theService->propagator("PropagatorWithMaterial")));
00426     result = fromInside(iniState);
00427   } else {
00428     StateOnTrackerBound fromOutside(&*propagator());
00429     result = fromOutside(innerState(staMuon));
00430   }
00431   return result;
00432 }
00433 
00434 TrajectorySeed FastTSGFromPropagation::createSeed(const TrajectoryStateOnSurface& tsos, const DetId& id) const {
00435 
00436   edm::OwnVector<TrackingRecHit> container;
00437   return createSeed(tsos, container, id);
00438 
00439 }
00440 
00441 TrajectorySeed FastTSGFromPropagation::createSeed(const TrajectoryStateOnSurface& tsos, const edm::OwnVector<TrackingRecHit>& container, const DetId& id) const {
00442 
00443   PTrajectoryStateOnDet seedTSOS = trajectoryStateTransform::persistentState(tsos,id.rawId());
00444   return TrajectorySeed(seedTSOS,container,oppositeToMomentum);
00445 
00446 }
00447 
00448 
00449 void FastTSGFromPropagation::validMeasurements(std::vector<TrajectoryMeasurement>& tms) const {
00450 
00451   std::vector<TrajectoryMeasurement>::iterator tmsend = std::remove_if(tms.begin(), tms.end(), isInvalid());
00452   tms.erase(tmsend, tms.end());
00453   return;
00454 
00455 }
00456 
00457 std::vector<TrajectoryMeasurement> FastTSGFromPropagation::findMeasurements_new(const DetLayer* nl, const TrajectoryStateOnSurface& staState) const {
00458 
00459   std::vector<TrajectoryMeasurement> result;
00460 
00461   std::vector<DetLayer::DetWithState> compatDets = nl->compatibleDets(staState, *propagator(), *estimator());
00462   if ( compatDets.empty() )  return result;
00463 
00464   for (std::vector<DetLayer::DetWithState>::const_iterator idws = compatDets.begin(); idws != compatDets.end(); ++idws) {
00465      if ( idws->second.isValid() && (idws->first) )  {
00466          std::vector<TrajectoryMeasurement> tmptm = 
00467            theMeasTracker->idToDet(idws->first->geographicalId())->fastMeasurements(idws->second, idws->second, *propagator(), *estimator());
00468          //validMeasurements(tmptm);
00469 //         if ( tmptm.size() > 2 ) {
00470 //            std::stable_sort(tmptm.begin(),tmptm.end(),increasingEstimate());
00471 //            result.insert(result.end(),tmptm.begin(), tmptm.begin()+2);
00472 //         } else {
00473             result.insert(result.end(),tmptm.begin(), tmptm.end());
00474 //         }
00475      }
00476   }
00477   
00478   return result;
00479 
00480 }
00481 
00482 std::vector<TrajectoryMeasurement> FastTSGFromPropagation::findMeasurements(const DetLayer* nl, const TrajectoryStateOnSurface& staState) const {
00483 
00484   std::vector<TrajectoryMeasurement> result = tkLayerMeasurements()->measurements((*nl), staState, *propagator(), *estimator());
00485   validMeasurements(result);
00486   return result;
00487 }
00488 
00489 bool FastTSGFromPropagation::passSelection(const TrajectoryStateOnSurface& tsos) const {
00490   if ( !theSelectStateFlag ) return true;
00491   else {
00492      if ( theBeamSpot.isValid() ) {
00493        return ( ( fabs(zDis(tsos) - theBeamSpot->z0() ) < theSigmaZ) );
00494 
00495      } else {
00496        return ( ( fabs(zDis(tsos)) < theSigmaZ) ); 
00497 //      double theDxyCut = 100;
00498 //      return ( (zDis(tsos) < theSigmaZ) && (dxyDis(tsos) < theDxyCut) );
00499      }
00500   }
00501 
00502 }
00503 
00504 double FastTSGFromPropagation::dxyDis(const TrajectoryStateOnSurface& tsos) const {
00505   return fabs(( - tsos.globalPosition().x() * tsos.globalMomentum().y() + tsos.globalPosition().y() * tsos.globalMomentum().x() )/tsos.globalMomentum().perp());
00506 }
00507 
00508 double FastTSGFromPropagation::zDis(const TrajectoryStateOnSurface& tsos) const {
00509   return tsos.globalPosition().z() - tsos.globalPosition().perp() * tsos.globalMomentum().z()/tsos.globalMomentum().perp();
00510 }
00511 
00512 void FastTSGFromPropagation::getRescalingFactor(const TrackCand& staMuon) {
00513     float pt = (staMuon.second)->pt();
00514     if ( pt < 13.0 ) theFlexErrorRescaling = 3; 
00515     else if ( pt < 30.0 ) theFlexErrorRescaling = 5;
00516     else theFlexErrorRescaling = 10;
00517     return;
00518 }
00519 
00520 
00521 void FastTSGFromPropagation::adjust(FreeTrajectoryState & state) const {
00522 
00523   //rescale the error
00524   if ( theResetMethod == "discreate" ) {
00525      state.rescaleError(theFlexErrorRescaling);
00526      return;
00527   }
00528 
00529   //rescale the error
00530   if ( theResetMethod == "fixed" || !theErrorMatrixAdjuster) {
00531      state.rescaleError(theFixedErrorRescaling);
00532      return;
00533   }
00534 
00535   CurvilinearTrajectoryError oMat = state.curvilinearError();
00536   CurvilinearTrajectoryError sfMat = theErrorMatrixAdjuster->get(state.momentum());//FIXME with position
00537   MuonErrorMatrix::multiply(oMat, sfMat);
00538   
00539   state = FreeTrajectoryState(state.parameters(),
00540                               oMat);
00541 }
00542 
00543 void FastTSGFromPropagation::adjust(TrajectoryStateOnSurface & state) const {
00544 
00545   //rescale the error
00546   if ( theResetMethod == "discreate" ) {
00547      state.rescaleError(theFlexErrorRescaling);
00548      return;
00549   }
00550 
00551   if ( theResetMethod == "fixed" || !theErrorMatrixAdjuster) {
00552      state.rescaleError(theFixedErrorRescaling);
00553      return;
00554   }
00555 
00556   CurvilinearTrajectoryError oMat = state.curvilinearError();
00557   CurvilinearTrajectoryError sfMat = theErrorMatrixAdjuster->get(state.globalMomentum());//FIXME with position
00558   MuonErrorMatrix::multiply(oMat, sfMat);
00559   
00560   state = TrajectoryStateOnSurface(state.globalParameters(),
00561                                    oMat,
00562                                    state.surface(),
00563                                    state.surfaceSide(),
00564                                    state.weight());
00565 }
00566 
00567 void FastTSGFromPropagation::stateOnDet(const TrajectoryStateOnSurface& ts,
00568                                         unsigned int detid,
00569                                         PTrajectoryStateOnDet& pts) const
00570 {
00571     const AlgebraicSymMatrix55& m = ts.localError().matrix();
00572     int dim = 5; 
00573     float localErrors[15];
00574     int k = 0;
00575     for (int i=0; i<dim; ++i) {
00576         for (int j=0; j<=i; ++j) {
00577             localErrors[k++] = m(i,j);
00578         }
00579     }
00580     int surfaceSide = static_cast<int>(ts.surfaceSide());
00581     pts = PTrajectoryStateOnDet( ts.localParameters(),
00582                                 localErrors, detid,
00583                                 surfaceSide);
00584 }