CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
FastTSGFromPropagation.cc
Go to the documentation of this file.
2 
18 
21 // #include "RecoTracker/MeasurementDet/interface/TkStripMeasurementDet.h"
25 
28 
32 
44 
45 using namespace std;
46 
47 
48 FastTSGFromPropagation::FastTSGFromPropagation(const edm::ParameterSet & iConfig,edm::ConsumesCollector& iC) :theTkLayerMeasurements(), theTracker(0), theNavigation(0), theService(0), theEstimator(0), theSigmaZ(0), theConfig (iConfig),
49  beamSpot_(iConfig.getParameter<edm::InputTag>("beamSpot"))
50 {
51  theCategory = "FastSimulation|Muons||FastTSGFromPropagation";
52  theMeasurementTrackerEventTag = iConfig.getParameter<edm::InputTag>("MeasurementTrackerEvent");
53 }
54 
55 FastTSGFromPropagation::FastTSGFromPropagation(const edm::ParameterSet & iConfig, const MuonServiceProxy* service,edm::ConsumesCollector& iC) : theTkLayerMeasurements(), theTracker(0), theNavigation(0), theService(service),theUpdator(0), theEstimator(0), theSigmaZ(0), theConfig (iConfig),
56  beamSpot_(iConfig.getParameter<edm::InputTag>("beamSpot"))
57 {
58  theCategory = "FastSimulation|Muons|FastTSGFromPropagation";
59  theMeasurementTrackerEventTag = iConfig.getParameter<edm::InputTag>("MeasurementTrackerEvent");
60 }
61 
63 {
64 
65  LogTrace(theCategory) << " FastTSGFromPropagation dtor called ";
66  if ( theNavigation ) delete theNavigation;
67  if ( theUpdator ) delete theUpdator;
68  if ( theEstimator ) delete theEstimator;
70 
71 }
72 
74  const TrackerTopology *tTopo, std::vector<TrajectorySeed> & result) {
75 
76  if ( theResetMethod == "discrete" ) getRescalingFactor(staMuon);
77 
78  TrajectoryStateOnSurface staState = outerTkState(staMuon);
79 
80  if ( !staState.isValid() ) {
81  LogTrace(theCategory) << "Error: initial state from L2 muon is invalid.";
82  return;
83  }
84 
85  LogTrace(theCategory) << "begin of trackerSeed:\n staState pos: "<<staState.globalPosition()
86  << " mom: "<<staState.globalMomentum()
87  <<"pos eta: "<<staState.globalPosition().eta()
88  <<"mom eta: "<<staState.globalMomentum().eta();
89 
90 
91  std::vector<const DetLayer*> nls = theNavigation->compatibleLayers(*(staState.freeState()), oppositeToMomentum);
92 
93  LogTrace(theCategory) << " compatible layers: "<<nls.size();
94 
95  if ( nls.empty() ) return;
96 
97  int ndesLayer = 0;
98 
99  bool usePredictedState = false;
100 
101  if ( theUpdateStateFlag ) { //use updated states
102  std::vector<TrajectoryMeasurement> alltm;
103 
104  for (std::vector<const DetLayer*>::const_iterator inl = nls.begin();
105  inl != nls.end(); inl++, ndesLayer++ ) {
106  if ( (*inl == 0) ) break;
107 // if ( (inl != nls.end()-1 ) && ( (*inl)->subDetector() == GeomDetEnumerators::TEC ) && ( (*(inl+1))->subDetector() == GeomDetEnumerators::TOB ) ) continue;
108  alltm = findMeasurements_new(*inl, staState);
109  if ( (!alltm.empty()) ) {
110  LogTrace(theCategory) << "final compatible layer: "<<ndesLayer;
111  break;
112  }
113  }
114 
115  if ( alltm.empty() ) {
116  LogTrace(theCategory) << " NO Measurements Found: eta: "<<staState.globalPosition().eta() <<"pt "<<staState.globalMomentum().perp();
117  usePredictedState = true;
118  } else {
119  LogTrace(theCategory) << " Measurements for seeds: "<<alltm.size();
120  std::stable_sort(alltm.begin(),alltm.end(),increasingEstimate());
121  if ( alltm.size() > 5 ) alltm.erase(alltm.begin() + 5, alltm.end());
122 
123  const edm::SimTrackContainer* simTracks = &(*theSimTracks);
124  const std::vector<unsigned> theSimTrackIds = theGSRecHits->ids();
125  TrackerRecHit theSeedHits;
126  std::vector<TrackerRecHit> outerHits;
127 
128  //std::vector<TrajectorySeed> tmpTS;
129  bool isMatch = false;
130  for (std::vector<TrajectoryMeasurement>::const_iterator itm = alltm.begin(); itm != alltm.end(); itm++) {
131  const TrajectoryStateOnSurface seedState = itm->predictedState();
132  double preY = seedState.globalPosition().y();
133 
134  // Check SimTrack
135  TrackingRecHit* aTrackingRecHit;
136  FreeTrajectoryState simtrack_trackerstate;
137  for( unsigned tkId=0; tkId != theSimTrackIds.size(); ++tkId ) {
138  const SimTrack & simtrack = (*simTracks)[theSimTrackIds[tkId]];
139  SiTrackerGSMatchedRecHit2DCollection::range theRecHitRange = theGSRecHits->get(theSimTrackIds[tkId]);
140  SiTrackerGSMatchedRecHit2DCollection::const_iterator theRecHitRangeIteratorBegin = theRecHitRange.first;
141  SiTrackerGSMatchedRecHit2DCollection::const_iterator theRecHitRangeIteratorEnd = theRecHitRange.second;
143 
145  simtrack.trackerSurfacePosition().y(),
146  simtrack.trackerSurfacePosition().z());
147  GlobalVector momentum(simtrack.trackerSurfaceMomentum().x(),
148  simtrack.trackerSurfaceMomentum().y(),
149  simtrack.trackerSurfaceMomentum().z());
150  int charge = (int)simtrack.charge();
151  GlobalTrajectoryParameters glb_parameters(position, momentum, charge, &*theService->magneticField().product());
152  simtrack_trackerstate = FreeTrajectoryState(glb_parameters);
153 
154  unsigned int outerId = 0;
155  for( iterRecHit = theRecHitRangeIteratorBegin; iterRecHit != theRecHitRangeIteratorEnd; ++iterRecHit) {
156  theSeedHits = TrackerRecHit(&(*iterRecHit), theGeometry, tTopo);
157  unsigned int id = theSeedHits.hit()->geographicalId().rawId();
158  if( preY < 0 ) {
159  if( id > outerId ) outerId = id;
160  }
161  else {
162  if( id > outerId ) outerId = id;
163  }
164  }
165  for( iterRecHit = theRecHitRangeIteratorBegin; iterRecHit != theRecHitRangeIteratorEnd; ++iterRecHit) {
166  theSeedHits = TrackerRecHit(&(*iterRecHit), theGeometry, tTopo);
167  if( itm->recHit()->hit()->geographicalId().rawId() == theSeedHits.hit()->geographicalId().rawId() ) {
168  aTrackingRecHit = theSeedHits.hit()->clone();
169  TransientTrackingRecHit::ConstRecHitPointer recHit = theTTRHBuilder->build(aTrackingRecHit);
170  if( !recHit ) continue;
171  TrajectoryStateOnSurface updatedTSOS = updator()->update(seedState, *(recHit));
172  if( updatedTSOS.isValid() && passSelection(updatedTSOS) ) {
174  container.push_back(recHit->hit()->clone());
175  TrajectorySeed ts = createSeed(updatedTSOS, container, recHit->geographicalId());
176  // check direction
177  const BasicTrajectorySeed* aSeed = &ts;
178  PTrajectoryStateOnDet PTSOD = aSeed->startingState();
179 
180  const GeomDet *g = theGeometry->idToDet(PTSOD.detId());
181  TrajectoryStateOnSurface tsos = trajectoryStateTransform::transientState(PTSOD, &(g->surface()), &*theService->magneticField().product());
182  if( tsos.globalMomentum().basicVector()*seedState.globalMomentum().basicVector() < 0. ) continue;
183  result.push_back(ts);
184  isMatch = true;
185  }
186  }
187  }
188  }
189  }
190  if( !isMatch ) {
191  // if there is no hits w.r.t. TM, find outermost hit
192  for (std::vector<TrajectoryMeasurement>::const_iterator itm = alltm.begin(); itm != alltm.end(); itm++) {
193  const TrajectoryStateOnSurface seedState = itm->predictedState();
194  double preY = seedState.globalPosition().y();
195 
196  // Check SimTrack
197  TrackingRecHit* aTrackingRecHit;
198  FreeTrajectoryState simtrack_trackerstate;
199  for( unsigned tkId=0; tkId != theSimTrackIds.size(); ++tkId ) {
200  const SimTrack & simtrack = (*simTracks)[theSimTrackIds[tkId]];
201  SiTrackerGSMatchedRecHit2DCollection::range theRecHitRange = theGSRecHits->get(theSimTrackIds[tkId]);
202  SiTrackerGSMatchedRecHit2DCollection::const_iterator theRecHitRangeIteratorBegin = theRecHitRange.first;
203  SiTrackerGSMatchedRecHit2DCollection::const_iterator theRecHitRangeIteratorEnd = theRecHitRange.second;
205 
207  simtrack.trackerSurfacePosition().y(),
208  simtrack.trackerSurfacePosition().z());
209  GlobalVector momentum(simtrack.trackerSurfaceMomentum().x(),
210  simtrack.trackerSurfaceMomentum().y(),
211  simtrack.trackerSurfaceMomentum().z());
212  int charge = (int)simtrack.charge();
213  GlobalTrajectoryParameters glb_parameters(position, momentum, charge, &*theService->magneticField().product());
214  simtrack_trackerstate = FreeTrajectoryState(glb_parameters);
215 
216  unsigned int outerId = 0;
217  for( iterRecHit = theRecHitRangeIteratorBegin; iterRecHit != theRecHitRangeIteratorEnd; ++iterRecHit) {
218  theSeedHits = TrackerRecHit(&(*iterRecHit), theGeometry, tTopo);
219  unsigned int id = theSeedHits.hit()->geographicalId().rawId();
220  if( preY < 0 ) {
221  if( id > outerId ) outerId = id;
222  }
223  else {
224  if( id > outerId ) outerId = id;
225  }
226  }
227  for( iterRecHit = theRecHitRangeIteratorBegin; iterRecHit != theRecHitRangeIteratorEnd; ++iterRecHit) {
228  theSeedHits = TrackerRecHit(&(*iterRecHit), theGeometry, tTopo);
229  if( outerId == theSeedHits.hit()->geographicalId().rawId() ) {
230  aTrackingRecHit = theSeedHits.hit()->clone();
231  TransientTrackingRecHit::ConstRecHitPointer recHit = theTTRHBuilder->build(aTrackingRecHit);
232  if( !recHit ) continue;
233  TrajectoryStateOnSurface updatedTSOS = updator()->update(seedState, *(recHit));
234  if( updatedTSOS.isValid() && passSelection(updatedTSOS) ) {
236  container.push_back(recHit->hit()->clone());
237  TrajectorySeed ts = createSeed(updatedTSOS, container, recHit->geographicalId());
238  // check direction
239  const BasicTrajectorySeed* aSeed = &ts;
240  PTrajectoryStateOnDet PTSOD = aSeed->startingState();
241 
242  const GeomDet *g = theGeometry->idToDet(PTSOD.detId());
243  TrajectoryStateOnSurface tsos = trajectoryStateTransform::transientState(PTSOD, &(g->surface()), &*theService->magneticField().product());
244  if( tsos.globalMomentum().basicVector()*seedState.globalMomentum().basicVector() < 0. ) continue;
245  result.push_back(ts);
246  }
247  }
248  }
249  }
250  }
251  }
252 
253  /*
254  for( unsigned ir = 0; ir < tmpTS.size(); ir++ ) {
255  const BasicTrajectorySeed* aSeed = &((tmpTS)[ir]);
256  PTrajectoryStateOnDet PTSOD = aSeed->startingState();
257 
258  DetId seedDetId(PTSOD.detId());
259  const GeomDet * g = theGeometry->idToDet(seedDetId);
260  TrajectoryStateOnSurface tsos = trajectoryStateTransform::transientState(PTSOD, &(g->surface()), &*theService->magneticField().product());
261  cout << "tsos3 = " << tsos.globalMomentum() << endl;
262  if( _index == ir ) {
263  cout << "tsos4 = " << tsos.globalMomentum() << endl;
264  result.push_back(tmpTS[ir]);
265  }
266  }
267  */
268  LogTrace(theCategory) << "result: "<<result.size();
269  return;
270  }
271  }
272 
273  if ( !theUpdateStateFlag || usePredictedState ) { //use predicted states
274  LogTrace(theCategory) << "use predicted state: ";
275  for (std::vector<const DetLayer*>::const_iterator inl = nls.begin();
276  inl != nls.end(); inl++ ) {
277 
278  if ( !result.empty() || *inl == 0 ) {
279  break;
280  }
281  std::vector<DetLayer::DetWithState> compatDets = (*inl)->compatibleDets(staState, *propagator(), *estimator());
282  LogTrace(theCategory) << " compatDets "<<compatDets.size();
283  if ( compatDets.empty() ) continue;
284  TrajectorySeed ts = createSeed(compatDets.front().second, compatDets.front().first->geographicalId());
285  result.push_back(ts);
286 
287  }
288  LogTrace(theCategory) << "result: "<<result.size();
289  return;
290  }
291  return;
292 }
293 
295 
296  theMaxChi2 = theConfig.getParameter<double>("MaxChi2");
297 
298  theFixedErrorRescaling = theConfig.getParameter<double>("ErrorRescaling");
299 
300  theFlexErrorRescaling = 1.0;
301 
303 
304  if (theResetMethod != "discrete" && theResetMethod != "fixed" && theResetMethod != "matrix" ) {
305  edm::LogError("FastTSGFromPropagation")
306  <<"Wrong error rescaling method: "<<theResetMethod <<"\n"
307  <<"Possible choices are: discrete, fixed, matrix.\n"
308  <<"Use discrete method" <<std::endl;
309  theResetMethod = "discrete";
310  }
311 
313 
314  theCacheId_MT = 0;
315 
316  theCacheId_TG = 0;
317 
319 
320  theService = service;
321 
322  theUseVertexStateFlag = theConfig.getParameter<bool>("UseVertexState");
323 
324  theUpdateStateFlag = theConfig.getParameter<bool>("UpdateState");
325 
326  theSelectStateFlag = theConfig.getParameter<bool>("SelectState");
327 
328  theSimTrackCollectionLabel = theConfig.getParameter<edm::InputTag>("SimTrackCollectionLabel");
330 
331  theUpdator = new KFUpdator();
332 
333  theSigmaZ = theConfig.getParameter<double>("SigmaZ");
334 
335  edm::ParameterSet errorMatrixPset = theConfig.getParameter<edm::ParameterSet>("errorMatrixPset");
336  if ( theResetMethod == "matrix" && !errorMatrixPset.empty()){
337  theAdjustAtIp = errorMatrixPset.getParameter<bool>("atIP");
338  theErrorMatrixAdjuster = new MuonErrorMatrix(errorMatrixPset);
339  } else {
340  theAdjustAtIp =false;
342  }
343 
344  theService->eventSetup().get<TrackerRecoGeometryRecord>().get(theTracker);
346 
348  theService->eventSetup().get<TrackerDigiGeometryRecord>().get(geometry);
349  theGeometry = &(*geometry);
350 
351  theService->eventSetup().get<TransientRecHitRecord>().get("WithTrackAngle", theTTRHBuilder);
352 
353 }
354 
356 
358 
359  // retrieve the MC truth (SimTracks)
362 
363 
364  unsigned long long newCacheId_MT = theService->eventSetup().get<CkfComponentsRecord>().cacheIdentifier();
365 
366  if ( theUpdateStateFlag && newCacheId_MT != theCacheId_MT ) {
367  LogTrace(theCategory) << "Measurment Tracker Geometry changed!";
368  theCacheId_MT = newCacheId_MT;
369  theService->eventSetup().get<CkfComponentsRecord>().get(theMeasTracker);
370  }
371 
372  if ( theUpdateStateFlag ) {
375  }
376 
377  bool trackerGeomChanged = false;
378 
379  unsigned long long newCacheId_TG = theService->eventSetup().get<TrackerRecoGeometryRecord>().cacheIdentifier();
380 
381  if ( newCacheId_TG != theCacheId_TG ) {
382  LogTrace(theCategory) << "Tracker Reco Geometry changed!";
383  theCacheId_TG = newCacheId_TG;
384  theService->eventSetup().get<TrackerRecoGeometryRecord>().get(theTracker);
385  trackerGeomChanged = true;
386  }
387 
388  if ( trackerGeomChanged && (&*theTracker) ) {
389  if ( theNavigation ) delete theNavigation;
391  }
392 }
393 
395 
396  TrajectoryStateOnSurface innerTS;
397 
398  if ( staMuon.first && staMuon.first->isValid() ) {
399  if (staMuon.first->direction() == alongMomentum) {
400  innerTS = staMuon.first->firstMeasurement().updatedState();
401  }
402  else if (staMuon.first->direction() == oppositeToMomentum) {
403  innerTS = staMuon.first->lastMeasurement().updatedState();
404  }
405  } else {
406  innerTS = trajectoryStateTransform::innerStateOnSurface(*(staMuon.second),*theService->trackingGeometry(), &*theService->magneticField());
407  }
408  //rescale the error
409  adjust(innerTS);
410 
411  return innerTS;
412 
413 }
414 
416 
418 
419  if ( theUseVertexStateFlag && staMuon.second->pt() > 1.0 ) {
420  FreeTrajectoryState iniState = trajectoryStateTransform::initialFreeState(*(staMuon.second), &*theService->magneticField());
421  //rescale the error at IP
422  adjust(iniState);
423 
424  StateOnTrackerBound fromInside(&*(theService->propagator("PropagatorWithMaterial")));
425  result = fromInside(iniState);
426  } else {
427  StateOnTrackerBound fromOutside(&*propagator());
428  result = fromOutside(innerState(staMuon));
429  }
430  return result;
431 }
432 
434 
436  return createSeed(tsos, container, id);
437 
438 }
439 
441 
443  return TrajectorySeed(seedTSOS,container,oppositeToMomentum);
444 
445 }
446 
447 
448 void FastTSGFromPropagation::validMeasurements(std::vector<TrajectoryMeasurement>& tms) const {
449 
450  std::vector<TrajectoryMeasurement>::iterator tmsend = std::remove_if(tms.begin(), tms.end(), isInvalid());
451  tms.erase(tmsend, tms.end());
452  return;
453 
454 }
455 
456 std::vector<TrajectoryMeasurement> FastTSGFromPropagation::findMeasurements_new(const DetLayer* nl, const TrajectoryStateOnSurface& staState) const {
457 
458  std::vector<TrajectoryMeasurement> result;
459 
460  std::vector<DetLayer::DetWithState> compatDets = nl->compatibleDets(staState, *propagator(), *estimator());
461  if ( compatDets.empty() ) return result;
462 
463  for (std::vector<DetLayer::DetWithState>::const_iterator idws = compatDets.begin(); idws != compatDets.end(); ++idws) {
464  if ( idws->second.isValid() && (idws->first) ) {
465  std::vector<TrajectoryMeasurement> tmptm =
466  theMeasTrackerEvent->idToDet(idws->first->geographicalId()).fastMeasurements(idws->second, idws->second, *propagator(), *estimator());
467  //validMeasurements(tmptm);
468 // if ( tmptm.size() > 2 ) {
469 // std::stable_sort(tmptm.begin(),tmptm.end(),increasingEstimate());
470 // result.insert(result.end(),tmptm.begin(), tmptm.begin()+2);
471 // } else {
472  result.insert(result.end(),tmptm.begin(), tmptm.end());
473 // }
474  }
475  }
476 
477  return result;
478 
479 }
480 
481 std::vector<TrajectoryMeasurement> FastTSGFromPropagation::findMeasurements(const DetLayer* nl, const TrajectoryStateOnSurface& staState) const {
482 
483  std::vector<TrajectoryMeasurement> result = tkLayerMeasurements()->measurements((*nl), staState, *propagator(), *estimator());
484  validMeasurements(result);
485  return result;
486 }
487 
489  if ( !theSelectStateFlag ) return true;
490  else {
491  if ( theBeamSpot.isValid() ) {
492  return ( ( fabs(zDis(tsos) - theBeamSpot->z0() ) < theSigmaZ) );
493 
494  } else {
495  return ( ( fabs(zDis(tsos)) < theSigmaZ) );
496 // double theDxyCut = 100;
497 // return ( (zDis(tsos) < theSigmaZ) && (dxyDis(tsos) < theDxyCut) );
498  }
499  }
500 
501 }
502 
504  return fabs(( - tsos.globalPosition().x() * tsos.globalMomentum().y() + tsos.globalPosition().y() * tsos.globalMomentum().x() )/tsos.globalMomentum().perp());
505 }
506 
508  return tsos.globalPosition().z() - tsos.globalPosition().perp() * tsos.globalMomentum().z()/tsos.globalMomentum().perp();
509 }
510 
512  float pt = (staMuon.second)->pt();
513  if ( pt < 13.0 ) theFlexErrorRescaling = 3;
514  else if ( pt < 30.0 ) theFlexErrorRescaling = 5;
515  else theFlexErrorRescaling = 10;
516  return;
517 }
518 
519 
521 
522  //rescale the error
523  if ( theResetMethod == "discreate" ) {
525  return;
526  }
527 
528  //rescale the error
529  if ( theResetMethod == "fixed" || !theErrorMatrixAdjuster) {
531  return;
532  }
533 
535  CurvilinearTrajectoryError sfMat = theErrorMatrixAdjuster->get(state.momentum());//FIXME with position
536  MuonErrorMatrix::multiply(oMat, sfMat);
537 
538  state = FreeTrajectoryState(state.parameters(),
539  oMat);
540 }
541 
543 
544  //rescale the error
545  if ( theResetMethod == "discreate" ) {
547  return;
548  }
549 
550  if ( theResetMethod == "fixed" || !theErrorMatrixAdjuster) {
552  return;
553  }
554 
556  CurvilinearTrajectoryError sfMat = theErrorMatrixAdjuster->get(state.globalMomentum());//FIXME with position
557  MuonErrorMatrix::multiply(oMat, sfMat);
558 
559  state = TrajectoryStateOnSurface(state.weight(),
560  state.globalParameters(),
561  oMat,
562  state.surface(),
563  state.surfaceSide());
564 }
565 
567  unsigned int detid,
568  PTrajectoryStateOnDet& pts) const
569 {
570  const AlgebraicSymMatrix55& m = ts.localError().matrix();
571  int dim = 5;
572  float localErrors[15];
573  int k = 0;
574  for (int i=0; i<dim; ++i) {
575  for (int j=0; j<=i; ++j) {
576  localErrors[k++] = m(i,j);
577  }
578  }
579  int surfaceSide = static_cast<int>(ts.surfaceSide());
581  localErrors, detid,
582  surfaceSide);
583 }
std::vector< TrajectoryMeasurement > findMeasurements(const DetLayer *, const TrajectoryStateOnSurface &) const
look for measurements on the first compatible layer
edm::Handle< MeasurementTrackerEvent > theMeasTrackerEvent
edm::ESHandle< MeasurementTracker > theMeasTracker
const math::XYZVectorD & trackerSurfacePosition() const
Definition: SimTrack.h:36
T getParameter(std::string const &) const
std::vector< TrajectoryMeasurement > measurements(const DetLayer &layer, const TrajectoryStateOnSurface &startingState, const Propagator &prop, const MeasurementEstimator &est) const
int i
Definition: DBlmapReader.cc:9
FastTSGFromPropagation(const edm::ParameterSet &pset, edm::ConsumesCollector &iC)
constructor
edm::Handle< reco::BeamSpot > theBeamSpot
std::pair< const_iterator, const_iterator > range
iterator range
Definition: RangeMap.h:50
T perp() const
Definition: PV3DBase.h:72
const Chi2MeasurementEstimator * theEstimator
void adjust(FreeTrajectoryState &) const
adjust the error matrix of the FTS
unsigned long long theCacheId_MT
std::pair< const Trajectory *, reco::TrackRef > TrackCand
TrajectoryStateOnSurface innerStateOnSurface(const reco::Track &tk, const TrackingGeometry &geom, const MagneticField *field)
const LocalTrajectoryParameters & localParameters() const
const LayerMeasurements * tkLayerMeasurements() const
const GlobalTrajectoryParameters & parameters() const
edm::Handle< edm::SimTrackContainer > theSimTracks
const DirectTrackerNavigation * theNavigation
const CurvilinearTrajectoryError & curvilinearError() const
virtual ~FastTSGFromPropagation()
destructor
edm::ESHandle< GeometricSearchTracker > theTracker
PTrajectoryStateOnDet persistentState(const TrajectoryStateOnSurface &ts, unsigned int detid)
T y() const
Definition: PV3DBase.h:63
GlobalPoint globalPosition() const
void validMeasurements(std::vector< TrajectoryMeasurement > &) const
select valid measurements
CurvilinearTrajectoryError get(GlobalVector momentum, bool convolute=true)
main method to be used. Retrieve a 5x5 symetrical matrix according to parametrization of error or sca...
const MuonServiceProxy * theService
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
float charge() const
charge
Definition: CoreSimTrack.cc:3
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:35
double charge(const std::vector< uint8_t > &Ampls)
const CurvilinearTrajectoryError & curvilinearError() const
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
virtual std::vector< DetWithState > compatibleDets(const TrajectoryStateOnSurface &startingState, const Propagator &prop, const MeasurementEstimator &est) const
static int position[TOTALCHAMBERS][3]
Definition: ReadPGInfo.cc:509
const Chi2MeasurementEstimator * estimator() const
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
void getRescalingFactor(const TrackCand &staMuon)
static const double pts[33]
Definition: Constants.h:30
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:43
void push_back(D *&d)
Definition: OwnVector.h:274
virtual TrajectoryStateOnSurface update(const TrajectoryStateOnSurface &, const TrackingRecHit &) const =0
int iEvent
Definition: GenABIO.cc:230
const SurfaceType & surface() const
edm::InputTag theMeasurementTrackerEventTag
bool passSelection(const TrajectoryStateOnSurface &) const
check some quantity and beam-spot compatibility and decide to continue
const TrajectoryStateUpdator * theUpdator
std::shared_ptr< TrackingRecHit const > ConstRecHitPointer
SurfaceSide surfaceSide() const
Position relative to material, defined relative to momentum vector.
FreeTrajectoryState const * freeState(bool withErrors=true) const
T z() const
Definition: PV3DBase.h:64
tuple result
Definition: query.py:137
int j
Definition: DBlmapReader.cc:9
std::vector< TrajectoryMeasurement > findMeasurements_new(const DetLayer *, const TrajectoryStateOnSurface &) const
look for measurements on the first compatible layer (faster way)
unsigned int detId() const
TrajectoryStateOnSurface outerTkState(const TrackCand &) const
const AlgebraicSymMatrix55 & matrix() const
static void multiply(CurvilinearTrajectoryError &initial_error, const CurvilinearTrajectoryError &scale_error)
multiply term by term the two matrix
TrajectoryStateOnSurface innerState(const TrackCand &) const
LayerMeasurements theTkLayerMeasurements
const LocalTrajectoryError & localError() const
bool isValid() const
Definition: HandleBase.h:76
virtual const GeomDet * idToDet(DetId) const
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:390
GlobalVector momentum() const
#define LogTrace(id)
int k[5][pyjets_maxn]
MuonErrorMatrix * theErrorMatrixAdjuster
edm::Handle< SiTrackerGSMatchedRecHit2DCollection > theGSRecHits
Definition: DetId.h:18
void rescaleError(double factor)
const TrackerGeometry * theGeometry
void init(const MuonServiceProxy *)
initialize
double zDis(const TrajectoryStateOnSurface &tsos) const
TrajectoryStateOnSurface transientState(const PTrajectoryStateOnDet &ts, const Surface *surface, const MagneticField *field)
const TrajectoryStateUpdator * updator() const
const GlobalTrajectoryParameters & globalParameters() const
FreeTrajectoryState initialFreeState(const reco::Track &tk, const MagneticField *field)
const math::XYZTLorentzVectorD & trackerSurfaceMomentum() const
Definition: SimTrack.h:38
edm::ESHandle< Propagator > propagator() const
TrajectorySeed createSeed(const TrajectoryStateOnSurface &, const DetId &) const
create a hitless seed from a trajectory state
void stateOnDet(const TrajectoryStateOnSurface &ts, unsigned int detid, PTrajectoryStateOnDet &pts) const
A mere copy (without memory leak) of an existing tracking method.
virtual GSSiTrackerRecHit2DLocalPos * clone() const =0
T eta() const
Definition: PV3DBase.h:76
const GSSiTrackerRecHit2DLocalPos * hit() const
Definition: TrackerRecHit.h:76
edm::InputTag theSimTrackCollectionLabel
ESHandle< TrackerGeometry > geometry
GlobalVector globalMomentum() const
std::vector< const DetLayer * > compatibleLayers(const FreeTrajectoryState &fts, PropagationDirection timeDirection) const
find compatible layers for a given trajectory state
double dxyDis(const TrajectoryStateOnSurface &tsos) const
DetId geographicalId() const
void trackerSeeds(const TrackCand &, const TrackingRegion &, const TrackerTopology *tTopo, std::vector< TrajectorySeed > &)
generate seed(s) for a track
T x() const
Definition: PV3DBase.h:62
const BasicVectorType & basicVector() const
Definition: PV3DBase.h:56
edm::ESHandle< TransientTrackingRecHitBuilder > theTTRHBuilder
std::vector< SimTrack > SimTrackContainer
void setEvent(const edm::Event &)
set an event
unsigned long long theCacheId_TG