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 NuclearInteractionFinder::NuclearInteractionFinder(const edm::EventSetup& es, const edm::ParameterSet& iConfig) :
00018 maxHits(iConfig.getParameter<int>("maxHits")),
00019 rescaleErrorFactor(iConfig.getParameter<double>("rescaleErrorFactor")),
00020 checkCompletedTrack(iConfig.getParameter<bool>("checkCompletedTrack")),
00021 navigationSchoolName(iConfig.getParameter<std::string>("NavigationSchool"))
00022 {
00023
00024 std::string measurementTrackerName = iConfig.getParameter<std::string>("MeasurementTrackerName");
00025
00026 edm::ESHandle<Propagator> prop;
00027 edm::ESHandle<TrajectoryStateUpdator> upd;
00028 edm::ESHandle<Chi2MeasurementEstimatorBase> est;
00029 edm::ESHandle<MeasurementTracker> measurementTrackerHandle;
00030 edm::ESHandle<GeometricSearchTracker> theGeomSearchTrackerHandle;
00031 edm::ESHandle<TrackerGeometry> theTrackerGeom;
00032
00033 es.get<TrackerDigiGeometryRecord> ().get (theTrackerGeom);
00034 es.get<TrackingComponentsRecord>().get("PropagatorWithMaterial",prop);
00035 es.get<TrackingComponentsRecord>().get("Chi2",est);
00036 es.get<CkfComponentsRecord>().get(measurementTrackerName, measurementTrackerHandle);
00037 es.get<TrackerRecoGeometryRecord>().get(theGeomSearchTrackerHandle );
00038 es.get<IdealMagneticFieldRecord>().get(theMagField);
00039
00040 edm::ESHandle<NavigationSchool> nav;
00041 es.get<NavigationSchoolRecord>().get(navigationSchoolName, nav);
00042 theNavigationSchool = nav.product();
00043
00044 thePropagator = prop.product();
00045 theEstimator = est.product();
00046 theMeasurementTracker = measurementTrackerHandle.product();
00047 theLayerMeasurements = new LayerMeasurements(theMeasurementTracker);
00048 theGeomSearchTracker = theGeomSearchTrackerHandle.product();
00049
00050 LogDebug("NuclearSeedGenerator") << "New NuclearInteractionFinder instance with parameters : \n"
00051 << "maxHits : " << maxHits << "\n"
00052 << "rescaleErrorFactor : " << rescaleErrorFactor << "\n"
00053 << "checkCompletedTrack : " << checkCompletedTrack << "\n";
00054
00055 nuclTester = new NuclearTester(maxHits, theEstimator, theTrackerGeom.product() );
00056
00057 currentSeed = new SeedFromNuclearInteraction(thePropagator, theTrackerGeom.product(), iConfig) ;
00058
00059 thePrimaryHelix = new TangentHelix();
00060 }
00061
00062 void NuclearInteractionFinder::setEvent(const edm::Event& event) const
00063 {
00064 theMeasurementTracker->update(event);
00065 }
00066
00067
00068 NuclearInteractionFinder::~NuclearInteractionFinder() {
00069 delete theLayerMeasurements;
00070 delete nuclTester;
00071 delete currentSeed;
00072 delete thePrimaryHelix;
00073 }
00074
00075
00076 bool NuclearInteractionFinder::run(const Trajectory& traj) {
00077
00078 if(traj.empty() || !traj.isValid()) return false;
00079
00080 LogDebug("NuclearSeedGenerator") << "Analyzis of a new trajectory with a number of valid hits = " << traj.foundHits();
00081
00082 std::vector<TrajectoryMeasurement> measurements = traj.measurements();
00083
00084
00085 nuclTester->reset( measurements.size() );
00086 allSeeds.clear();
00087
00088
00089 if(traj.direction()==alongMomentum) {
00090 std::reverse(measurements.begin(), measurements.end());
00091 }
00092
00093 std::vector<TrajectoryMeasurement>::const_iterator it_meas = measurements.begin();
00094
00095 std::vector<double> ncompatibleHits;
00096 bool NIfound = false;
00097
00098
00099 while(!NIfound)
00100 {
00101 if(it_meas == measurements.end()) break;
00102
00103
00104 if(nuclTester->nHitsChecked() > maxHits) break;
00105
00106 nuclTester->push_back(*it_meas, findCompatibleMeasurements(*it_meas, rescaleErrorFactor));
00107
00108 LogDebug("NuclearSeedGenerator") << "Number of compatible meas:" << (nuclTester->back()).size() << "\n"
00109 << "Mean distance between hits :" << nuclTester->meanHitDistance() << "\n"
00110 << "Forward estimate :" << nuclTester->fwdEstimate() << "\n";
00111
00112
00113 if( checkCompletedTrack==false && (nuclTester->compatibleHits()).front()==0 ) break;
00114
00115 if(nuclTester->isNuclearInteraction()) NIfound=true;
00116
00117 ++it_meas;
00118 }
00119
00120 if(NIfound) {
00121 LogDebug("NuclearSeedGenerator") << "NUCLEAR INTERACTION FOUND at index : " << nuclTester->nuclearIndex() << "\n";
00122
00123
00124 definePrimaryHelix(measurements.begin()+nuclTester->nuclearIndex()-1);
00125
00126 this->fillSeeds( nuclTester->goodTMPair());
00127
00128 return true;
00129 }
00130
00131 return false;
00132 }
00133
00134 void NuclearInteractionFinder::definePrimaryHelix(std::vector<TrajectoryMeasurement>::const_iterator it_meas) {
00135
00136
00137 GlobalPoint pt[3];
00138 for(int i=0; i<3; i++) {
00139 pt[i] = (it_meas->updatedState()).globalParameters().position();
00140 it_meas++;
00141 }
00142 delete thePrimaryHelix;
00143 thePrimaryHelix = new TangentHelix( pt[0], pt[1], pt[2] );
00144 }
00145
00146 std::vector<TrajectoryMeasurement>
00147 NuclearInteractionFinder::findCompatibleMeasurements(const TM& lastMeas, double rescale) const
00148 {
00149 TSOS currentState = lastMeas.updatedState();
00150 LogDebug("NuclearSeedGenerator") << "currentState :" << currentState << "\n";
00151
00152 TSOS newState = rescaleError(rescale, currentState);
00153 return findMeasurementsFromTSOS(newState, lastMeas.recHit()->geographicalId());
00154 }
00155
00156 std::vector<TrajectoryMeasurement>
00157 NuclearInteractionFinder::findMeasurementsFromTSOS(const TSOS& currentState, DetId detid) const {
00158
00159 using namespace std;
00160 int invalidHits = 0;
00161 vector<TM> result;
00162 const DetLayer* lastLayer = theGeomSearchTracker->detLayer( detid );
00163 vector<const DetLayer*> nl;
00164
00165 if(lastLayer) {
00166 nl = lastLayer->nextLayers( *currentState.freeState(), alongMomentum);
00167 }
00168 else {
00169 edm::LogError("NuclearInteractionFinder") << "In findCompatibleMeasurements : lastLayer not accessible";
00170 return result;
00171 }
00172
00173 if (nl.empty()) {
00174 LogDebug("NuclearSeedGenerator") << "In findCompatibleMeasurements : no compatible layer found";
00175 return result;
00176 }
00177
00178 for (vector<const DetLayer*>::iterator il = nl.begin();
00179 il != nl.end(); il++) {
00180 vector<TM> tmp =
00181 theLayerMeasurements->measurements((**il),currentState, *thePropagator, *theEstimator);
00182 if ( !tmp.empty()) {
00183 if ( result.empty()) result = tmp;
00184 else {
00185
00186 result.insert( result.end()-invalidHits, tmp.begin(), tmp.end());
00187 }
00188 invalidHits++;
00189 }
00190 }
00191
00192
00193 if ( result.size() > 1) {
00194 sort( result.begin(), result.end()-invalidHits, TrajMeasLessEstim());
00195 }
00196 return result;
00197 }
00198
00199
00200 void NuclearInteractionFinder::fillSeeds( const std::pair<TrajectoryMeasurement, std::vector<TrajectoryMeasurement> >& tmPairs ) {
00201
00202
00203 const TM& innerTM = tmPairs.first;
00204 const std::vector<TM>& outerTMs = tmPairs.second;
00205
00206
00207 for(std::vector<TM>::const_iterator outtm = outerTMs.begin(); outtm!=outerTMs.end(); outtm++) {
00208 if((innerTM.recHit())->isValid() && (outtm->recHit())->isValid()) {
00209 currentSeed->setMeasurements(innerTM.updatedState(), innerTM.recHit(), outtm->recHit());
00210 allSeeds.push_back(*currentSeed);
00211 }
00212 else LogDebug("NuclearSeedGenerator") << "The initial hits for seeding are invalid" << "\n";
00213 }
00214 return;
00215 }
00216
00217 std::auto_ptr<TrajectorySeedCollection> NuclearInteractionFinder::getPersistentSeeds() {
00218 std::auto_ptr<TrajectorySeedCollection> output(new TrajectorySeedCollection);
00219 for(std::vector<SeedFromNuclearInteraction>::const_iterator it_seed = allSeeds.begin(); it_seed != allSeeds.end(); it_seed++) {
00220 if(it_seed->isValid()) {
00221 output->push_back( it_seed->TrajSeed() );
00222 }
00223 else LogDebug("NuclearSeedGenerator") << "The seed is invalid" << "\n";
00224 }
00225 return output;
00226 }
00227
00228 void NuclearInteractionFinder::improveSeeds() {
00229 std::vector<SeedFromNuclearInteraction> newSeedCollection;
00230
00231
00232 for(std::vector<SeedFromNuclearInteraction>::const_iterator it_seed = allSeeds.begin(); it_seed != allSeeds.end(); it_seed++) {
00233
00234 if( !it_seed->isValid() ) continue;
00235
00236
00237 std::vector<TM> thirdTMs = findMeasurementsFromTSOS( it_seed->updatedTSOS() , it_seed->outerHitDetId() );
00238
00239
00240 for(std::vector<TM>::const_iterator tm = thirdTMs.begin(); tm!= thirdTMs.end(); tm++) {
00241
00242 if( ! tm->recHit()->isValid() ) continue;
00243
00244
00245 currentSeed->setMeasurements(*thePrimaryHelix, it_seed->initialTSOS(), it_seed->outerHit(), tm->recHit() );
00246 newSeedCollection.push_back( *currentSeed );
00247 }
00248 }
00249 allSeeds.clear();
00250 allSeeds = newSeedCollection;
00251 }
00252
00253 TrajectoryStateOnSurface NuclearInteractionFinder::rescaleError(float rescale, const TSOS& state) const {
00254
00255 AlgebraicSymMatrix55 m(state.localError().matrix());
00256 AlgebraicSymMatrix55 mr;
00257 LocalTrajectoryParameters ltp = state.localParameters();
00258
00259
00260 mr(0,0) = (ltp.signedInverseMomentum()*0.2*rescale)*(ltp.signedInverseMomentum()*0.2*rescale);
00261
00262
00263 mr(1,1) = 1E-2*rescale*rescale;
00264 mr(2,2) = 1E-2*rescale*rescale;
00265
00266
00267 mr(3,3) = m(3,3);
00268 mr(4,4) = m(4,4);
00269
00270 return TSOS(ltp, mr, state.surface(), &(state.globalParameters().magneticField()), state.surfaceSide());
00271 }