CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_4/src/RecoTracker/NuclearSeedGenerator/src/NuclearInteractionFinder.cc

Go to the documentation of this file.
00001 #include "RecoTracker/NuclearSeedGenerator/interface/NuclearInteractionFinder.h"
00002 #include "TrackingTools/TrackFitters/interface/TrajectoryStateWithArbitraryError.h"
00003 #include "TrackingTools/PatternTools/interface/TrajMeasLessEstim.h"
00004 
00005 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00006 
00007 #include "FWCore/Utilities/interface/Exception.h"
00008 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00009 
00010 #include "TrackingTools/KalmanUpdators/interface/Chi2MeasurementEstimatorBase.h"
00011 #include "RecoTracker/Record/interface/NavigationSchoolRecord.h" 
00012 
00013 #include "RecoTracker/Record/interface/CkfComponentsRecord.h"
00014 
00015 #include "DataFormats/TrajectorySeed/interface/PropagationDirection.h"
00016 
00017 #include "FWCore/Framework/interface/ESHandle.h"
00018 
00019 
00020 NuclearInteractionFinder::NuclearInteractionFinder(const edm::EventSetup& es, const edm::ParameterSet& iConfig) :
00021 maxHits(iConfig.getParameter<int>("maxHits")),
00022 rescaleErrorFactor(iConfig.getParameter<double>("rescaleErrorFactor")),
00023 checkCompletedTrack(iConfig.getParameter<bool>("checkCompletedTrack")),
00024 navigationSchoolName(iConfig.getParameter<std::string>("NavigationSchool"))
00025 {
00026 
00027    std::string measurementTrackerName = iConfig.getParameter<std::string>("MeasurementTrackerName");
00028 
00029    edm::ESHandle<Propagator> prop;
00030    edm::ESHandle<TrajectoryStateUpdator> upd;
00031    edm::ESHandle<Chi2MeasurementEstimatorBase> est;
00032    edm::ESHandle<MeasurementTracker> measurementTrackerHandle;
00033    edm::ESHandle<GeometricSearchTracker>       theGeomSearchTrackerHandle;
00034    edm::ESHandle<TrackerGeometry>              theTrackerGeom;
00035 
00036    es.get<TrackerDigiGeometryRecord> ().get (theTrackerGeom);
00037    es.get<TrackingComponentsRecord>().get("PropagatorWithMaterial",prop);
00038    es.get<TrackingComponentsRecord>().get("Chi2",est);
00039    es.get<CkfComponentsRecord>().get(measurementTrackerName, measurementTrackerHandle);
00040    es.get<TrackerRecoGeometryRecord>().get(theGeomSearchTrackerHandle );
00041    es.get<IdealMagneticFieldRecord>().get(theMagField);
00042 
00043    edm::ESHandle<NavigationSchool>             nav;
00044    es.get<NavigationSchoolRecord>().get(navigationSchoolName, nav);
00045    theNavigationSchool  = nav.product();
00046 
00047    thePropagator = prop.product();
00048    theEstimator = est.product();
00049    theMeasurementTracker = measurementTrackerHandle.product();
00050    theLayerMeasurements = new LayerMeasurements(theMeasurementTracker);
00051    theGeomSearchTracker = theGeomSearchTrackerHandle.product();
00052 
00053    LogDebug("NuclearSeedGenerator") << "New NuclearInteractionFinder instance with parameters : \n"
00054                                         << "maxHits : " << maxHits << "\n"
00055                                         << "rescaleErrorFactor : " << rescaleErrorFactor << "\n"
00056                                         << "checkCompletedTrack : " << checkCompletedTrack << "\n";
00057 
00058    nuclTester = new NuclearTester(maxHits, theEstimator, theTrackerGeom.product() );
00059 
00060    currentSeed = new SeedFromNuclearInteraction(thePropagator, theTrackerGeom.product(), iConfig) ;
00061 
00062    thePrimaryHelix = new TangentHelix();
00063 }
00064 //----------------------------------------------------------------------
00065 void NuclearInteractionFinder::setEvent(const edm::Event& event) const
00066 {
00067    theMeasurementTracker->update(event);
00068 }
00069 
00070 //----------------------------------------------------------------------
00071 NuclearInteractionFinder::~NuclearInteractionFinder() {
00072   delete theLayerMeasurements;
00073   delete nuclTester;
00074   delete currentSeed;
00075   delete thePrimaryHelix;
00076 }
00077 
00078 //----------------------------------------------------------------------
00079 bool  NuclearInteractionFinder::run(const Trajectory& traj) {
00080 
00081         if(traj.empty() || !traj.isValid()) return false;
00082 
00083         LogDebug("NuclearSeedGenerator") << "Analyzis of a new trajectory with a number of valid hits = " << traj.foundHits();
00084 
00085         std::vector<TrajectoryMeasurement> measurements = traj.measurements();
00086 
00087         // initialization
00088         nuclTester->reset( measurements.size() );
00089         allSeeds.clear();
00090 
00091 
00092         if(traj.direction()==alongMomentum)  {
00093                 std::reverse(measurements.begin(), measurements.end());
00094         }
00095 
00096         std::vector<TrajectoryMeasurement>::const_iterator it_meas = measurements.begin();
00097 
00098         std::vector<double> ncompatibleHits;
00099         bool NIfound = false;
00100 
00101         // Loop on all the RecHits.
00102         while(!NIfound)
00103          {
00104            if(it_meas == measurements.end()) break;
00105 
00106            // check only the maxHits outermost hits of the primary track
00107            if(nuclTester->nHitsChecked() > maxHits) break;
00108 
00109            nuclTester->push_back(*it_meas, findCompatibleMeasurements(*it_meas, rescaleErrorFactor));
00110 
00111            LogDebug("NuclearSeedGenerator") << "Number of compatible meas:" << (nuclTester->back()).size() << "\n"
00112                                                 << "Mean distance between hits :" << nuclTester->meanHitDistance() << "\n"
00113                                                 << "Forward estimate :" << nuclTester->fwdEstimate() << "\n";
00114 
00115            // don't check tracks which reach the end of the tracker if checkCompletedTrack==false
00116            if( checkCompletedTrack==false && (nuclTester->compatibleHits()).front()==0 ) break;
00117 
00118            if(nuclTester->isNuclearInteraction()) NIfound=true;
00119 
00120            ++it_meas;
00121          }
00122 
00123         if(NIfound) {
00124             LogDebug("NuclearSeedGenerator") << "NUCLEAR INTERACTION FOUND at index : " << nuclTester->nuclearIndex()  << "\n";
00125 
00126             // Get correct parametrization of the helix of the primary track at the interaction point (to be used by improveCurrentSeed)
00127             definePrimaryHelix(measurements.begin()+nuclTester->nuclearIndex()-1);
00128 
00129             this->fillSeeds( nuclTester->goodTMPair());
00130 
00131             return true;
00132         }
00133 
00134     return false;
00135 }
00136 //----------------------------------------------------------------------
00137 void NuclearInteractionFinder::definePrimaryHelix(std::vector<TrajectoryMeasurement>::const_iterator it_meas) {
00138     // This method uses the 3 last TM after the interaction point to calculate the helix parameters
00139 
00140     GlobalPoint pt[3];
00141     for(int i=0; i<3; i++) {
00142        pt[i] = (it_meas->updatedState()).globalParameters().position();
00143        it_meas++;
00144     }
00145     delete thePrimaryHelix;
00146     thePrimaryHelix = new TangentHelix( pt[0], pt[1], pt[2] );
00147 }
00148 //----------------------------------------------------------------------
00149 std::vector<TrajectoryMeasurement>
00150 NuclearInteractionFinder::findCompatibleMeasurements(const TM& lastMeas, double rescale) const
00151 {
00152   TSOS currentState = lastMeas.updatedState();
00153   LogDebug("NuclearSeedGenerator") << "currentState :" << currentState << "\n";
00154 
00155   TSOS newState = rescaleError(rescale, currentState);
00156   return findMeasurementsFromTSOS(newState, lastMeas.recHit()->geographicalId());
00157 }
00158 //----------------------------------------------------------------------
00159 std::vector<TrajectoryMeasurement>
00160 NuclearInteractionFinder::findMeasurementsFromTSOS(const TSOS& currentState, DetId detid) const {
00161 
00162   using namespace std;
00163   int invalidHits = 0;
00164   vector<TM> result;
00165   const DetLayer* lastLayer = theGeomSearchTracker->detLayer( detid );
00166   vector<const DetLayer*> nl;
00167 
00168   if(lastLayer) {
00169           nl = lastLayer->nextLayers( *currentState.freeState(), alongMomentum);
00170   }
00171   else {
00172       edm::LogError("NuclearInteractionFinder") << "In findCompatibleMeasurements : lastLayer not accessible";
00173       return result;
00174   }
00175 
00176   if (nl.empty()) {
00177       LogDebug("NuclearSeedGenerator") << "In findCompatibleMeasurements :  no compatible layer found";
00178       return result;
00179   }
00180 
00181   for (vector<const DetLayer*>::iterator il = nl.begin();
00182        il != nl.end(); il++) {
00183     vector<TM> tmp =
00184       theLayerMeasurements->measurements((**il),currentState, *thePropagator, *theEstimator);
00185     if ( !tmp.empty()) {
00186       if ( result.empty()) result = tmp;
00187       else {
00188         // keep one dummy TM at the end, skip the others
00189         result.insert( result.end()-invalidHits, tmp.begin(), tmp.end());
00190       }
00191       invalidHits++;
00192     }
00193   }
00194 
00195   // sort the final result, keep dummy measurements at the end
00196   if ( result.size() > 1) {
00197     sort( result.begin(), result.end()-invalidHits, TrajMeasLessEstim());
00198   }
00199   return result;
00200 }
00201 
00202 //----------------------------------------------------------------------
00203 void NuclearInteractionFinder::fillSeeds( const std::pair<TrajectoryMeasurement, std::vector<TrajectoryMeasurement> >& tmPairs ) {
00204   // This method returns the seeds calculated by the class SeedsFromNuclearInteraction
00205 
00206             const TM& innerTM = tmPairs.first;
00207             const std::vector<TM>& outerTMs = tmPairs.second;
00208 
00209             // Loop on all outer TM
00210             for(std::vector<TM>::const_iterator outtm = outerTMs.begin(); outtm!=outerTMs.end(); outtm++) {
00211                if((innerTM.recHit())->isValid() && (outtm->recHit())->isValid()) {
00212                      currentSeed->setMeasurements(innerTM.updatedState(), innerTM.recHit(), outtm->recHit());
00213                      allSeeds.push_back(*currentSeed);
00214                 }
00215                 else  LogDebug("NuclearSeedGenerator") << "The initial hits for seeding are invalid" << "\n";
00216              }
00217              return;
00218 }
00219 //----------------------------------------------------------------------
00220 std::auto_ptr<TrajectorySeedCollection> NuclearInteractionFinder::getPersistentSeeds() {
00221    std::auto_ptr<TrajectorySeedCollection> output(new TrajectorySeedCollection);
00222    for(std::vector<SeedFromNuclearInteraction>::const_iterator it_seed = allSeeds.begin(); it_seed != allSeeds.end(); it_seed++) {
00223        if(it_seed->isValid()) {
00224            output->push_back( it_seed->TrajSeed() );
00225        }
00226        else LogDebug("NuclearSeedGenerator") << "The seed is invalid" << "\n";
00227    }
00228    return output;
00229 }
00230 //----------------------------------------------------------------------
00231 void NuclearInteractionFinder::improveSeeds() {
00232         std::vector<SeedFromNuclearInteraction> newSeedCollection;
00233 
00234         // loop on all actual seeds
00235         for(std::vector<SeedFromNuclearInteraction>::const_iterator it_seed = allSeeds.begin(); it_seed != allSeeds.end(); it_seed++) {
00236 
00237               if( !it_seed->isValid() ) continue;
00238 
00239               // find compatible TM in an outer layer
00240               std::vector<TM> thirdTMs = findMeasurementsFromTSOS( it_seed->updatedTSOS() , it_seed->outerHitDetId() );
00241 
00242               // loop on those new TMs
00243               for(std::vector<TM>::const_iterator tm = thirdTMs.begin(); tm!= thirdTMs.end(); tm++) {
00244 
00245                    if( ! tm->recHit()->isValid() ) continue;
00246 
00247                    // create new seeds collection using the circle equation
00248                    currentSeed->setMeasurements(*thePrimaryHelix, it_seed->initialTSOS(), it_seed->outerHit(), tm->recHit() );
00249                    newSeedCollection.push_back( *currentSeed );
00250               }
00251        }
00252        allSeeds.clear();
00253        allSeeds = newSeedCollection;
00254 }
00255 //----------------------------------------------------------------------
00256 TrajectoryStateOnSurface NuclearInteractionFinder::rescaleError(float rescale, const TSOS& state) const {
00257 
00258      AlgebraicSymMatrix55 m(state.localError().matrix());
00259      AlgebraicSymMatrix55 mr;
00260      LocalTrajectoryParameters ltp = state.localParameters();
00261 
00262      // we assume that the error on q/p is equal to 20% of q/p * rescale 
00263      mr(0,0) = (ltp.signedInverseMomentum()*0.2*rescale)*(ltp.signedInverseMomentum()*0.2*rescale);
00264 
00265      // the error on dx/z and dy/dz is fixed to 10% (* rescale)
00266      mr(1,1) = 1E-2*rescale*rescale;
00267      mr(2,2) = 1E-2*rescale*rescale;
00268 
00269      // the error on the local x and y positions are not modified.
00270      mr(3,3) = m(3,3);
00271      mr(4,4) = m(4,4);
00272 
00273      return TSOS(ltp, mr, state.surface(), &(state.globalParameters().magneticField()), state.surfaceSide());
00274 }