![]() |
![]() |
#include <RecoTracker/NuclearSeedGenerator/interface/NuclearInteractionFinder.h>
Public Member Functions | |
std::auto_ptr < TrajectorySeedCollection > | getPersistentSeeds () |
Fill 'output' with persistent nuclear seeds. | |
void | improveSeeds () |
Improve the seeds with a third RecHit. | |
const NavigationSchool * | nav () const |
NuclearInteractionFinder (const edm::EventSetup &es, const edm::ParameterSet &iConfig) | |
NuclearInteractionFinder () | |
TrajectoryStateOnSurface | rescaleError (float rescale, const TSOS &state) const |
bool | run (const Trajectory &traj) |
Run the Finder. | |
void | setEvent (const edm::Event &event) const |
define the measurement (to be called for each event) | |
virtual | ~NuclearInteractionFinder () |
Private Types | |
typedef TrajectoryMeasurement::ConstRecHitPointer | ConstRecHitPointer |
typedef FreeTrajectoryState | FTS |
typedef TrajectoryMeasurement | TM |
typedef std::vector< Trajectory > | TrajectoryContainer |
typedef TrajectoryStateOnSurface | TSOS |
Private Member Functions | |
void | definePrimaryHelix (std::vector< TrajectoryMeasurement >::const_iterator it_meas) |
Calculate the parameters of the circle representing the primary track at the interaction point. | |
void | fillSeeds (const std::pair< TrajectoryMeasurement, std::vector< TrajectoryMeasurement > > &tmPairs) |
get the seeds at the interaction point | |
std::vector < TrajectoryMeasurement > | findCompatibleMeasurements (const TM &lastMeas, double rescaleFactor) const |
Find compatible TM of a TM with error rescaled by rescaleFactor. | |
std::vector < TrajectoryMeasurement > | findMeasurementsFromTSOS (const TSOS ¤tState, DetId detid) const |
Private Attributes | |
std::vector < SeedFromNuclearInteraction > | allSeeds |
bool | checkCompletedTrack |
If set to true check all the tracks, even those reaching the edge of the tracker. | |
SeedFromNuclearInteraction * | currentSeed |
unsigned int | maxHits |
std::string | navigationSchoolName |
NuclearTester * | nuclTester |
double | ptMin |
double | rescaleErrorFactor |
const MeasurementEstimator * | theEstimator |
const GeometricSearchTracker * | theGeomSearchTracker |
const LayerMeasurements * | theLayerMeasurements |
edm::ESHandle< MagneticField > | theMagField |
const MeasurementTracker * | theMeasurementTracker |
const NavigationSchool * | theNavigationSchool |
TangentHelix * | thePrimaryHelix |
const Propagator * | thePropagator |
The method run gets all compatible TMs of all TMs associated of a trajectory. Then it uses the NuclearTester class to decide whether the trajectory has interacted nuclearly or not. It finally returns a pair of the TM where the nuclear interaction occurs and all compatible TMs associated.
Definition at line 44 of file NuclearInteractionFinder.h.
typedef TrajectoryMeasurement::ConstRecHitPointer NuclearInteractionFinder::ConstRecHitPointer [private] |
Definition at line 51 of file NuclearInteractionFinder.h.
typedef FreeTrajectoryState NuclearInteractionFinder::FTS [private] |
Definition at line 48 of file NuclearInteractionFinder.h.
typedef TrajectoryMeasurement NuclearInteractionFinder::TM [private] |
Definition at line 49 of file NuclearInteractionFinder.h.
typedef std::vector<Trajectory> NuclearInteractionFinder::TrajectoryContainer [private] |
Definition at line 50 of file NuclearInteractionFinder.h.
typedef TrajectoryStateOnSurface NuclearInteractionFinder::TSOS [private] |
Definition at line 47 of file NuclearInteractionFinder.h.
NuclearInteractionFinder::NuclearInteractionFinder | ( | ) | [inline] |
NuclearInteractionFinder::NuclearInteractionFinder | ( | const edm::EventSetup & | es, | |
const edm::ParameterSet & | iConfig | |||
) |
Definition at line 17 of file NuclearInteractionFinder.cc.
References checkCompletedTrack, currentSeed, edm::EventSetup::get(), edm::ParameterSet::getParameter(), LogDebug, maxHits, HLT_VtxMuL3::measurementTrackerName, nav(), navigationSchoolName, nuclTester, edm::ESHandle< T >::product(), rescaleErrorFactor, theEstimator, theGeomSearchTracker, theLayerMeasurements, theMagField, theMeasurementTracker, theNavigationSchool, thePrimaryHelix, and thePropagator.
00017 : 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 } //----------------------------------------------------------------------
NuclearInteractionFinder::~NuclearInteractionFinder | ( | ) | [virtual] |
Definition at line 68 of file NuclearInteractionFinder.cc.
References currentSeed, nuclTester, theLayerMeasurements, and thePrimaryHelix.
00068 { 00069 delete theLayerMeasurements; 00070 delete nuclTester; 00071 delete currentSeed; 00072 delete thePrimaryHelix; 00073 }
void NuclearInteractionFinder::definePrimaryHelix | ( | std::vector< TrajectoryMeasurement >::const_iterator | it_meas | ) | [private] |
Calculate the parameters of the circle representing the primary track at the interaction point.
Definition at line 134 of file NuclearInteractionFinder.cc.
References i, and thePrimaryHelix.
Referenced by run().
00134 { 00135 // This method uses the 3 last TM after the interaction point to calculate the helix parameters 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 }
void NuclearInteractionFinder::fillSeeds | ( | const std::pair< TrajectoryMeasurement, std::vector< TrajectoryMeasurement > > & | tmPairs | ) | [private] |
get the seeds at the interaction point
Definition at line 200 of file NuclearInteractionFinder.cc.
References allSeeds, currentSeed, LogDebug, TrajectoryMeasurement::recHit(), SeedFromNuclearInteraction::setMeasurements(), and TrajectoryMeasurement::updatedState().
Referenced by run().
00200 { 00201 // This method returns the seeds calculated by the class SeedsFromNuclearInteraction 00202 00203 const TM& innerTM = tmPairs.first; 00204 const std::vector<TM>& outerTMs = tmPairs.second; 00205 00206 // Loop on all outer TM 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 }
std::vector< TrajectoryMeasurement > NuclearInteractionFinder::findCompatibleMeasurements | ( | const TM & | lastMeas, | |
double | rescaleFactor | |||
) | const [private] |
Find compatible TM of a TM with error rescaled by rescaleFactor.
Definition at line 147 of file NuclearInteractionFinder.cc.
References findMeasurementsFromTSOS(), LogDebug, TrajectoryMeasurement::recHit(), rescaleError(), and TrajectoryMeasurement::updatedState().
Referenced by run().
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 }
std::vector< TrajectoryMeasurement > NuclearInteractionFinder::findMeasurementsFromTSOS | ( | const TSOS & | currentState, | |
DetId | detid | |||
) | const [private] |
Definition at line 157 of file NuclearInteractionFinder.cc.
References alongMomentum, GeometricSearchTracker::detLayer(), TrajectoryStateOnSurface::freeState(), LogDebug, LayerMeasurements::measurements(), DetLayer::nextLayers(), HLT_VtxMuL3::result, python::multivaluedict::sort(), std, theEstimator, theGeomSearchTracker, theLayerMeasurements, thePropagator, and tmp.
Referenced by findCompatibleMeasurements(), and improveSeeds().
00157 { 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 // keep one dummy TM at the end, skip the others 00186 result.insert( result.end()-invalidHits, tmp.begin(), tmp.end()); 00187 } 00188 invalidHits++; 00189 } 00190 } 00191 00192 // sort the final result, keep dummy measurements at the end 00193 if ( result.size() > 1) { 00194 sort( result.begin(), result.end()-invalidHits, TrajMeasLessEstim()); 00195 } 00196 return result; 00197 }
std::auto_ptr< TrajectorySeedCollection > NuclearInteractionFinder::getPersistentSeeds | ( | ) |
Fill 'output' with persistent nuclear seeds.
Definition at line 217 of file NuclearInteractionFinder.cc.
References allSeeds, LogDebug, and output().
00217 { 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 }
void NuclearInteractionFinder::improveSeeds | ( | ) |
Improve the seeds with a third RecHit.
Definition at line 228 of file NuclearInteractionFinder.cc.
References allSeeds, currentSeed, findMeasurementsFromTSOS(), SeedFromNuclearInteraction::setMeasurements(), and thePrimaryHelix.
00228 { 00229 std::vector<SeedFromNuclearInteraction> newSeedCollection; 00230 00231 // loop on all actual seeds 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 // find compatible TM in an outer layer 00237 std::vector<TM> thirdTMs = findMeasurementsFromTSOS( it_seed->updatedTSOS() , it_seed->outerHitDetId() ); 00238 00239 // loop on those new TMs 00240 for(std::vector<TM>::const_iterator tm = thirdTMs.begin(); tm!= thirdTMs.end(); tm++) { 00241 00242 if( ! tm->recHit()->isValid() ) continue; 00243 00244 // create new seeds collection using the circle equation 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 }
const NavigationSchool* NuclearInteractionFinder::nav | ( | ) | const [inline] |
Definition at line 88 of file NuclearInteractionFinder.h.
References theNavigationSchool.
Referenced by NuclearInteractionFinder().
00088 { return theNavigationSchool; }
TrajectoryStateOnSurface NuclearInteractionFinder::rescaleError | ( | float | rescale, | |
const TSOS & | state | |||
) | const |
Definition at line 253 of file NuclearInteractionFinder.cc.
References TrajectoryStateOnSurface::globalParameters(), TrajectoryStateOnSurface::localError(), TrajectoryStateOnSurface::localParameters(), m, GlobalTrajectoryParameters::magneticField(), LocalTrajectoryError::matrix(), LocalTrajectoryParameters::signedInverseMomentum(), TrajectoryStateOnSurface::surface(), and TrajectoryStateOnSurface::surfaceSide().
Referenced by findCompatibleMeasurements().
00253 { 00254 00255 AlgebraicSymMatrix55 m(state.localError().matrix()); 00256 AlgebraicSymMatrix55 mr; 00257 LocalTrajectoryParameters ltp = state.localParameters(); 00258 00259 // we assume that the error on q/p is equal to 20% of q/p * rescale 00260 mr(0,0) = (ltp.signedInverseMomentum()*0.2*rescale)*(ltp.signedInverseMomentum()*0.2*rescale); 00261 00262 // the error on dx/z and dy/dz is fixed to 10% (* rescale) 00263 mr(1,1) = 1E-2*rescale*rescale; 00264 mr(2,2) = 1E-2*rescale*rescale; 00265 00266 // the error on the local x and y positions are not modified. 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 }
bool NuclearInteractionFinder::run | ( | const Trajectory & | traj | ) |
Run the Finder.
Definition at line 76 of file NuclearInteractionFinder.cc.
References allSeeds, alongMomentum, NuclearTester::back(), checkCompletedTrack, NuclearTester::compatibleHits(), definePrimaryHelix(), Trajectory::direction(), Trajectory::empty(), fillSeeds(), findCompatibleMeasurements(), Trajectory::foundHits(), prof2calltree::front, NuclearTester::fwdEstimate(), NuclearTester::goodTMPair(), NuclearTester::isNuclearInteraction(), Trajectory::isValid(), LogDebug, maxHits, NuclearTester::meanHitDistance(), Trajectory::measurements(), NuclearTester::nHitsChecked(), NuclearTester::nuclearIndex(), nuclTester, NuclearTester::push_back(), rescaleErrorFactor, NuclearTester::reset(), and size.
00076 { 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 // initialization 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 // Loop on all the RecHits. 00099 while(!NIfound) 00100 { 00101 if(it_meas == measurements.end()) break; 00102 00103 // check only the maxHits outermost hits of the primary track 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 // don't check tracks which reach the end of the tracker if checkCompletedTrack==false 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 // Get correct parametrization of the helix of the primary track at the interaction point (to be used by improveCurrentSeed) 00124 definePrimaryHelix(measurements.begin()+nuclTester->nuclearIndex()-1); 00125 00126 this->fillSeeds( nuclTester->goodTMPair()); 00127 00128 return true; 00129 } 00130 00131 return false; 00132 }
void NuclearInteractionFinder::setEvent | ( | const edm::Event & | event | ) | const |
define the measurement (to be called for each event)
Definition at line 62 of file NuclearInteractionFinder.cc.
References theMeasurementTracker, and MeasurementTracker::update().
00063 { 00064 theMeasurementTracker->update(event); 00065 }
std::vector< SeedFromNuclearInteraction > NuclearInteractionFinder::allSeeds [private] |
Definition at line 102 of file NuclearInteractionFinder.h.
Referenced by fillSeeds(), getPersistentSeeds(), improveSeeds(), and run().
If set to true check all the tracks, even those reaching the edge of the tracker.
Definition at line 109 of file NuclearInteractionFinder.h.
Referenced by NuclearInteractionFinder(), and run().
Definition at line 101 of file NuclearInteractionFinder.h.
Referenced by fillSeeds(), improveSeeds(), NuclearInteractionFinder(), and ~NuclearInteractionFinder().
unsigned int NuclearInteractionFinder::maxHits [private] |
Definition at line 107 of file NuclearInteractionFinder.h.
Referenced by NuclearInteractionFinder(), and run().
std::string NuclearInteractionFinder::navigationSchoolName [private] |
Definition at line 110 of file NuclearInteractionFinder.h.
Referenced by NuclearInteractionFinder().
Definition at line 100 of file NuclearInteractionFinder.h.
Referenced by NuclearInteractionFinder(), run(), and ~NuclearInteractionFinder().
double NuclearInteractionFinder::ptMin [private] |
Definition at line 106 of file NuclearInteractionFinder.h.
double NuclearInteractionFinder::rescaleErrorFactor [private] |
Definition at line 108 of file NuclearInteractionFinder.h.
Referenced by NuclearInteractionFinder(), and run().
const MeasurementEstimator* NuclearInteractionFinder::theEstimator [private] |
Definition at line 93 of file NuclearInteractionFinder.h.
Referenced by findMeasurementsFromTSOS(), and NuclearInteractionFinder().
const GeometricSearchTracker* NuclearInteractionFinder::theGeomSearchTracker [private] |
Definition at line 96 of file NuclearInteractionFinder.h.
Referenced by findMeasurementsFromTSOS(), and NuclearInteractionFinder().
const LayerMeasurements* NuclearInteractionFinder::theLayerMeasurements [private] |
Definition at line 95 of file NuclearInteractionFinder.h.
Referenced by findMeasurementsFromTSOS(), NuclearInteractionFinder(), and ~NuclearInteractionFinder().
const MeasurementTracker* NuclearInteractionFinder::theMeasurementTracker [private] |
Definition at line 94 of file NuclearInteractionFinder.h.
Referenced by NuclearInteractionFinder(), and setEvent().
const NavigationSchool* NuclearInteractionFinder::theNavigationSchool [private] |
Definition at line 97 of file NuclearInteractionFinder.h.
Referenced by nav(), and NuclearInteractionFinder().
Definition at line 103 of file NuclearInteractionFinder.h.
Referenced by definePrimaryHelix(), improveSeeds(), NuclearInteractionFinder(), and ~NuclearInteractionFinder().
const Propagator* NuclearInteractionFinder::thePropagator [private] |
Definition at line 92 of file NuclearInteractionFinder.h.
Referenced by findMeasurementsFromTSOS(), and NuclearInteractionFinder().