test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
NuclearInteractionFinder.cc
Go to the documentation of this file.
4 
6 
9 
12 
14 
16 
18 
19 
21 maxHits(iConfig.getParameter<int>("maxHits")),
22 rescaleErrorFactor(iConfig.getParameter<double>("rescaleErrorFactor")),
23 checkCompletedTrack(iConfig.getParameter<bool>("checkCompletedTrack")),
24 navigationSchoolName(iConfig.getParameter<std::string>("NavigationSchool"))
25 {
26 
27  std::string measurementTrackerName = iConfig.getParameter<std::string>("MeasurementTrackerName");
28 
32  edm::ESHandle<MeasurementTracker> measurementTrackerHandle;
33  edm::ESHandle<GeometricSearchTracker> theGeomSearchTrackerHandle;
34  edm::ESHandle<TrackerGeometry> theTrackerGeom;
35 
36  es.get<TrackerDigiGeometryRecord> ().get (theTrackerGeom);
37  es.get<TrackingComponentsRecord>().get("PropagatorWithMaterial",prop);
38  es.get<TrackingComponentsRecord>().get("Chi2",est);
39  es.get<CkfComponentsRecord>().get(measurementTrackerName, measurementTrackerHandle);
40  es.get<TrackerRecoGeometryRecord>().get(theGeomSearchTrackerHandle );
42 
46 
47  thePropagator = prop.product();
48  theEstimator = est.product();
49  theMeasurementTracker = measurementTrackerHandle.product();
50  theGeomSearchTracker = theGeomSearchTrackerHandle.product();
51 
52  LogDebug("NuclearSeedGenerator") << "New NuclearInteractionFinder instance with parameters : \n"
53  << "maxHits : " << maxHits << "\n"
54  << "rescaleErrorFactor : " << rescaleErrorFactor << "\n"
55  << "checkCompletedTrack : " << checkCompletedTrack << "\n";
56 
57  nuclTester = new NuclearTester(maxHits, theEstimator, theTrackerGeom.product() );
58 
59  currentSeed = new SeedFromNuclearInteraction(thePropagator, theTrackerGeom.product(), iConfig) ;
60 
62 }
63 //----------------------------------------------------------------------
65  delete nuclTester;
66  delete currentSeed;
67  delete thePrimaryHelix;
68 }
69 
70 //----------------------------------------------------------------------
72 
73  if(traj.empty() || !traj.isValid()) return false;
74 
75  LogDebug("NuclearSeedGenerator") << "Analyzis of a new trajectory with a number of valid hits = " << traj.foundHits();
76 
77  std::vector<TrajectoryMeasurement> measurements = traj.measurements();
78 
79  // initialization
80  nuclTester->reset( measurements.size() );
81  allSeeds.clear();
82 
83  LayerMeasurements layerMeasurements(*theMeasurementTracker, event);
84 
85  if(traj.direction()==alongMomentum) {
86  std::reverse(measurements.begin(), measurements.end());
87  }
88 
89  std::vector<TrajectoryMeasurement>::const_iterator it_meas = measurements.begin();
90 
91  std::vector<double> ncompatibleHits;
92  bool NIfound = false;
93 
94  // Loop on all the RecHits.
95  while(!NIfound)
96  {
97  if(it_meas == measurements.end()) break;
98 
99  // check only the maxHits outermost hits of the primary track
100  if(nuclTester->nHitsChecked() > maxHits) break;
101 
102  nuclTester->push_back(*it_meas, findCompatibleMeasurements(*it_meas, rescaleErrorFactor, layerMeasurements));
103 
104  LogDebug("NuclearSeedGenerator") << "Number of compatible meas:" << (nuclTester->back()).size() << "\n"
105  << "Mean distance between hits :" << nuclTester->meanHitDistance() << "\n"
106  << "Forward estimate :" << nuclTester->fwdEstimate() << "\n";
107 
108  // don't check tracks which reach the end of the tracker if checkCompletedTrack==false
109  if( checkCompletedTrack==false && (nuclTester->compatibleHits()).front()==0 ) break;
110 
111  if(nuclTester->isNuclearInteraction()) NIfound=true;
112 
113  ++it_meas;
114  }
115 
116  if(NIfound) {
117  LogDebug("NuclearSeedGenerator") << "NUCLEAR INTERACTION FOUND at index : " << nuclTester->nuclearIndex() << "\n";
118 
119  // Get correct parametrization of the helix of the primary track at the interaction point (to be used by improveCurrentSeed)
120  definePrimaryHelix(measurements.begin()+nuclTester->nuclearIndex()-1);
121 
122  this->fillSeeds( nuclTester->goodTMPair());
123 
124  return true;
125  }
126 
127  return false;
128 }
129 //----------------------------------------------------------------------
130 void NuclearInteractionFinder::definePrimaryHelix(std::vector<TrajectoryMeasurement>::const_iterator it_meas) {
131  // This method uses the 3 last TM after the interaction point to calculate the helix parameters
132 
133  GlobalPoint pt[3];
134  for(int i=0; i<3; i++) {
135  pt[i] = (it_meas->updatedState()).globalParameters().position();
136  it_meas++;
137  }
138  delete thePrimaryHelix;
139  thePrimaryHelix = new TangentHelix( pt[0], pt[1], pt[2] );
140 }
141 //----------------------------------------------------------------------
142 std::vector<TrajectoryMeasurement>
143 NuclearInteractionFinder::findCompatibleMeasurements(const TM& lastMeas, double rescale, const LayerMeasurements & layerMeasurements) const
144 {
145  TSOS currentState = lastMeas.updatedState();
146  LogDebug("NuclearSeedGenerator") << "currentState :" << currentState << "\n";
147 
148  TSOS newState = rescaleError(rescale, currentState);
149  return findMeasurementsFromTSOS(newState, lastMeas.recHit()->geographicalId(), layerMeasurements);
150 }
151 //----------------------------------------------------------------------
152 std::vector<TrajectoryMeasurement>
153 NuclearInteractionFinder::findMeasurementsFromTSOS(const TSOS& currentState, DetId detid, const LayerMeasurements & layerMeasurements) const {
154 
155  using namespace std;
156  int invalidHits = 0;
157  vector<TM> result;
158  const DetLayer* lastLayer = theGeomSearchTracker->detLayer( detid );
159  vector<const DetLayer*> nl;
160 
161  if(lastLayer) {
162  nl = theNavigationSchool->nextLayers(*lastLayer,*currentState.freeState(), alongMomentum);
163  }
164  else {
165  edm::LogError("NuclearInteractionFinder") << "In findCompatibleMeasurements : lastLayer not accessible";
166  return result;
167  }
168 
169  if (nl.empty()) {
170  LogDebug("NuclearSeedGenerator") << "In findCompatibleMeasurements : no compatible layer found";
171  return result;
172  }
173 
174  for (vector<const DetLayer*>::iterator il = nl.begin();
175  il != nl.end(); il++) {
176  vector<TM> tmp =
177  layerMeasurements.measurements((**il),currentState, *thePropagator, *theEstimator);
178  if ( !tmp.empty()) {
179  if ( result.empty()) result = tmp;
180  else {
181  // keep one dummy TM at the end, skip the others
182  result.insert( result.end()-invalidHits, tmp.begin(), tmp.end());
183  }
184  invalidHits++;
185  }
186  }
187 
188  // sort the final result, keep dummy measurements at the end
189  if ( result.size() > 1) {
190  sort( result.begin(), result.end()-invalidHits, TrajMeasLessEstim());
191  }
192  return result;
193 }
194 
195 //----------------------------------------------------------------------
196 void NuclearInteractionFinder::fillSeeds( const std::pair<TrajectoryMeasurement, std::vector<TrajectoryMeasurement> >& tmPairs ) {
197  // This method returns the seeds calculated by the class SeedsFromNuclearInteraction
198 
199  const TM& innerTM = tmPairs.first;
200  const std::vector<TM>& outerTMs = tmPairs.second;
201 
202  // Loop on all outer TM
203  for(std::vector<TM>::const_iterator outtm = outerTMs.begin(); outtm!=outerTMs.end(); outtm++) {
204  if((innerTM.recHit())->isValid() && (outtm->recHit())->isValid()) {
205  currentSeed->setMeasurements(innerTM.updatedState(), innerTM.recHit(), outtm->recHit());
206  allSeeds.push_back(*currentSeed);
207  }
208  else LogDebug("NuclearSeedGenerator") << "The initial hits for seeding are invalid" << "\n";
209  }
210  return;
211 }
212 //----------------------------------------------------------------------
213 std::auto_ptr<TrajectorySeedCollection> NuclearInteractionFinder::getPersistentSeeds() {
214  std::auto_ptr<TrajectorySeedCollection> output(new TrajectorySeedCollection);
215  for(std::vector<SeedFromNuclearInteraction>::const_iterator it_seed = allSeeds.begin(); it_seed != allSeeds.end(); it_seed++) {
216  if(it_seed->isValid()) {
217  output->push_back( it_seed->TrajSeed() );
218  }
219  else LogDebug("NuclearSeedGenerator") << "The seed is invalid" << "\n";
220  }
221  return output;
222 }
223 //----------------------------------------------------------------------
225  std::vector<SeedFromNuclearInteraction> newSeedCollection;
226 
227  LayerMeasurements layerMeasurements(*theMeasurementTracker, event);
228 
229  // loop on all actual seeds
230  for(std::vector<SeedFromNuclearInteraction>::const_iterator it_seed = allSeeds.begin(); it_seed != allSeeds.end(); it_seed++) {
231 
232  if( !it_seed->isValid() ) continue;
233 
234  // find compatible TM in an outer layer
235  std::vector<TM> thirdTMs = findMeasurementsFromTSOS( it_seed->updatedTSOS() , it_seed->outerHitDetId(), layerMeasurements );
236 
237  // loop on those new TMs
238  for(std::vector<TM>::const_iterator tm = thirdTMs.begin(); tm!= thirdTMs.end(); tm++) {
239 
240  if( ! tm->recHit()->isValid() ) continue;
241 
242  // create new seeds collection using the circle equation
243  currentSeed->setMeasurements(*thePrimaryHelix, it_seed->initialTSOS(), it_seed->outerHit(), tm->recHit() );
244  newSeedCollection.push_back( *currentSeed );
245  }
246  }
247  allSeeds.clear();
248  allSeeds = newSeedCollection;
249 }
250 //----------------------------------------------------------------------
252 
256 
257  // we assume that the error on q/p is equal to 20% of q/p * rescale
258  mr(0,0) = (ltp.signedInverseMomentum()*0.2*rescale)*(ltp.signedInverseMomentum()*0.2*rescale);
259 
260  // the error on dx/z and dy/dz is fixed to 10% (* rescale)
261  mr(1,1) = 1E-2*rescale*rescale;
262  mr(2,2) = 1E-2*rescale*rescale;
263 
264  // the error on the local x and y positions are not modified.
265  mr(3,3) = m(3,3);
266  mr(4,4) = m(4,4);
267 
268  return TSOS(ltp, mr, state.surface(), &(state.globalParameters().magneticField()), state.surfaceSide());
269 }
#define LogDebug(id)
bool empty() const
True if trajectory has no measurements.
Definition: Trajectory.h:299
T getParameter(std::string const &) const
int foundHits() const
Definition: Trajectory.h:279
std::vector< TrajectoryMeasurement > measurements(const DetLayer &layer, const TrajectoryStateOnSurface &startingState, const Propagator &prop, const MeasurementEstimator &est) const
int i
Definition: DBlmapReader.cc:9
ConstRecHitPointer const & recHit() const
std::vector< int > compatibleHits() const
Definition: NuclearTester.h:61
const LocalTrajectoryParameters & localParameters() const
const TMContainer & back() const
Definition: NuclearTester.h:43
bool run(const Trajectory &traj, const MeasurementTrackerEvent &event)
Run the Finder.
std::vector< TrajectoryMeasurement > findMeasurementsFromTSOS(const TSOS &currentState, DetId detid, const LayerMeasurements &layerMeasurements) const
void improveSeeds(const MeasurementTrackerEvent &event)
Improve the seeds with a third RecHit.
const NavigationSchool * nav() const
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
void setMeasurements(const TSOS &tsosAtInteractionPoint, ConstRecHitPointer ihit, ConstRecHitPointer ohit)
Fill all data members from 2 TM&#39;s where the first one is supposed to be at the interaction point...
double fwdEstimate(const std::vector< TrajectoryMeasurement > &vecTM) const
int nuclearIndex() const
Definition: NuclearTester.h:55
tuple result
Definition: mps_fire.py:84
PropagationDirection const & direction() const
Definition: Trajectory.cc:125
bool isNuclearInteraction()
float signedInverseMomentum() const
Signed inverse momentum q/p (zero for neutrals).
DataContainer const & measurements() const
Definition: Trajectory.h:250
void push_back(const TM &init_tm, const TMContainer &vecTM)
Definition: NuclearTester.h:38
const SurfaceType & surface() const
std::vector< TrajectorySeed > TrajectorySeedCollection
std::vector< TrajectoryMeasurement > findCompatibleMeasurements(const TM &lastMeas, double rescaleFactor, const LayerMeasurements &layerMeasurements) const
Find compatible TM of a TM with error rescaled by rescaleFactor.
SurfaceSide surfaceSide() const
Position relative to material, defined relative to momentum vector.
double meanHitDistance(const std::vector< TrajectoryMeasurement > &vecTM) const
FreeTrajectoryState const * freeState(bool withErrors=true) const
std::vector< SeedFromNuclearInteraction > allSeeds
void reset(unsigned int nMeasurements)
Definition: NuclearTester.h:49
const MeasurementTracker * theMeasurementTracker
const DetLayer * detLayer(const DetId &id) const
obsolete method. Use idToLayer() instead.
const NavigationSchool * theNavigationSchool
const AlgebraicSymMatrix55 & matrix() const
const MeasurementEstimator * theEstimator
const LocalTrajectoryError & localError() const
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
Definition: DetId.h:18
bool isValid() const
Definition: Trajectory.h:324
const GlobalTrajectoryParameters & globalParameters() const
std::vector< const DetLayer * > nextLayers(const DetLayer &detLayer, Args &&...args) const
NavigationDirection.
const T & get() const
Definition: EventSetup.h:56
edm::ESHandle< MagneticField > theMagField
T const * product() const
Definition: ESHandle.h:86
const GeometricSearchTracker * theGeomSearchTracker
TrajectoryStateOnSurface TSOS
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
std::auto_ptr< TrajectorySeedCollection > getPersistentSeeds()
Fill &#39;output&#39; with persistent nuclear seeds.
void fillSeeds(const std::pair< TrajectoryMeasurement, std::vector< TrajectoryMeasurement > > &tmPairs)
get the seeds at the interaction point
Class used to test if a track has interacted nuclearly.
Definition: NuclearTester.h:16
const MagneticField & magneticField() const
TrajectoryStateOnSurface const & updatedState() const
unsigned int nHitsChecked() const
Definition: NuclearTester.h:59
void definePrimaryHelix(std::vector< TrajectoryMeasurement >::const_iterator it_meas)
Calculate the parameters of the circle representing the primary track at the interaction point...
TrajectoryStateOnSurface rescaleError(float rescale, const TSOS &state) const
const TMPair & goodTMPair() const
Definition: NuclearTester.h:57
tuple size
Write out results.
SeedFromNuclearInteraction * currentSeed