CMS 3D CMS Logo

DynamicTruncation.cc
Go to the documentation of this file.
1 
30 
31 #define MAX_THR 1e7
32 
33 using namespace edm;
34 using namespace std;
35 using namespace reco;
36 
37 namespace dyt_utils{
38  static const std::map<etaRegion, std::string> etaRegionStr { {etaRegion::eta0p8, "eta0p8"},
39  {etaRegion::eta1p2, "eta1p2"},
40  {etaRegion::eta2p0, "eta2p0"},
41  {etaRegion::eta2p2, "eta2p2"},
42  {etaRegion::eta2p4, "eta2p4"}};
43 };
44 
45 
46 
48  propagator = theService.propagator("SmartPropagatorAny");
49  propagatorPF = theService.propagator("SmartPropagatorAny");
50  propagatorCompatibleDet = theService.propagator("SmartPropagatorAny");
51  theG = theService.trackingGeometry();
52  theService.eventSetup().get<TransientRecHitRecord>().get("MuonRecHitBuilder",theMuonRecHitBuilder);
53  theService.eventSetup().get<TrackingComponentsRecord>().get("KFUpdator",updatorHandle);
54  theService.eventSetup().get<MuonGeometryRecord>().get(cscGeom);
55  theService.eventSetup().get<MuonRecoGeometryRecord>().get(navMuon);
56  theService.eventSetup().get<IdealMagneticFieldRecord>().get(magfield);
57  navigation = new DirectMuonNavigation(theService.detLayerGeometry());
58  getSegs = new ChamberSegmentUtility();
59  thrManager = new ThrParameters(&theService.eventSetup());
60  useDBforThr = thrManager->isValidThdDB();
61  if (useDBforThr) dytThresholds = thrManager->getInitialThresholds();
62 
63  doUpdateOfKFStates = true;
64  useParametrizedThr = false;
65 }
66 
68  delete navigation;
69  delete thrManager;
70  delete getSegs;
71 }
72 
74  TrajectoryStateOnSurface temp = updatorHandle->update(tsos, *rechit);
75  if (temp.isValid()) tsos = updatorHandle->update(tsos, *rechit);
76 }
77 
79  ConstRecHitContainer tmprecHits;
80  vector<const TrackingRecHit*> DTrh = bestDTSeg.recHits();
81  for (vector<const TrackingRecHit*>::iterator it = DTrh.begin(); it != DTrh.end(); it++) {
82  tmprecHits.push_back(theMuonRecHitBuilder->build(*it));
83  }
84  sort(tmprecHits);
85  for (ConstRecHitContainer::const_iterator it = tmprecHits.begin(); it != tmprecHits.end(); ++it) {
86  DTLayerId layid((*it)->det()->geographicalId());
87  TrajectoryStateOnSurface temp = propagator->propagate(tsos, theG->idToDet(layid)->surface());
88  if (temp.isValid()) {
89  TrajectoryStateOnSurface tempTsos = updatorHandle->update(temp, **it);
90  if (tempTsos.isValid() ) tsos = tempTsos;
91  }
92  }
93 }
94 
96  ConstRecHitContainer tmprecHits;
97  vector<CSCRecHit2D> CSCrh = bestCSCSeg.specificRecHits();
98  for (vector<CSCRecHit2D>::iterator it = CSCrh.begin(); it != CSCrh.end(); ++it) {
99  tmprecHits.push_back(theMuonRecHitBuilder->build(&*it));
100  }
101  sort(tmprecHits);
102  for (ConstRecHitContainer::const_iterator it = tmprecHits.begin(); it != tmprecHits.end(); ++it) {
103  const CSCLayer* cscLayer = cscGeom->layer((*it)->det()->geographicalId());
104  TrajectoryStateOnSurface temp = propagator->propagate(tsos, cscLayer->surface());
105  if (temp.isValid()) {
106  TrajectoryStateOnSurface tempTsos = updatorHandle->update(temp, **it);
107  if (tempTsos.isValid() ) tsos = tempTsos;
108  }
109  }
110 }
111 
112 
116 void DynamicTruncation::setSelector(int selector) {
117  if (selector < 0 || selector > 2) throw cms::Exception("NotAvailable") << "DYT selector: wrong option!" << endl;
118  //if (selector == 0) cout << "[DYT disabled]\n";
119  //if (selector == 1) cout << "[use all compatible stations]\n";
120  //if (selector == 2) cout << "[stop at second consecutive incompatible station]\n";
121  DYTselector = selector;
122 
123 }
124 
125 void DynamicTruncation::setUseAPE(bool useAPE_) {
126  useAPE = useAPE_;
127 }
128 
130  doUpdateOfKFStates = upState;
131 }
132 
133 void DynamicTruncation::setThr(const vector<int>& thr) {
134  if (thr.size() == 2) {
135  for (unsigned int i = 0; i < thr.size(); i++)
136  if (thr[i] >= 0) Thrs.push_back(thr[i]);
137  else Thrs.push_back(MAX_THR);
138  return;
139  }
140  throw cms::Exception("NotAvailable") << "WARNING: wrong size for the threshold vector!\nExpected size: 2\n Found size: " << thr.size();
141 }
142 
144  for (auto const& region : dyt_utils::etaRegionStr ){
145  parameters[region.first] = par.getParameter< std::vector<double> >(region.second);
146  }
147 }
151 
152 
153 
154 //===> filter
156  result.clear();
157  prelFitMeas.clear();
158 
159  // Get APE maps
160  dtApeMap = thrManager->GetDTApeMap();
161  cscApeMap = thrManager->GetCSCApeMap();
162 
163  // Get Last tracker TSOS (updated)
164  vector<TrajectoryMeasurement> muonMeasurements = traj.measurements();
165  TrajectoryMeasurement lastTKm = muonMeasurements.front();
166  for (vector<TrajectoryMeasurement>::const_iterator imT = muonMeasurements.begin(); imT != muonMeasurements.end(); imT++ ) {
167  if ( !(*imT).recHit()->isValid() ) continue;
168  const TransientTrackingRecHit* hit = &(*(*imT).recHit());
169  if (hit->geographicalId().det() == DetId::Tracker) {
170  result.push_back((*imT).recHit());
171  if (!(*imT).forwardPredictedState().isValid()) continue;
172  if ((*imT).forwardPredictedState().globalPosition().mag() >
173  lastTKm.forwardPredictedState().globalPosition().mag()) lastTKm = *imT;
174  }
175  }
176  currentState = lastTKm.forwardPredictedState();
177  update(currentState, lastTKm.recHit());
178 
179  prelFitState = lastTKm.forwardPredictedState();
180  update(prelFitState, lastTKm.recHit());
181  prelFitMeas = result;
182 
183  // Run the DYT
184  filteringAlgo();
185 
186  return result;
187 }
188 
189 
190 //===> filteringAlgo
192  map<int, vector<DetId> > compatibleIds;
193  map<int, vector<DTRecSegment4D> > dtSegMap;
194  map<int, vector<CSCSegment> > cscSegMap;
195  int incompConLay = 0;
196  nStationsUsed = 0;
197 
198  // Get list of compatible layers
199  compatibleDets(currentState, compatibleIds);
200 
201  // Fill segment maps
202  fillSegmentMaps(compatibleIds, dtSegMap, cscSegMap);
203 
204  // Do a preliminary fit
205  if (useDBforThr) preliminaryFit(compatibleIds, dtSegMap, cscSegMap);
206 
207  // Loop on compatible layers
208  for (map<int, vector<DetId> >::iterator it=compatibleIds.begin(); it!=compatibleIds.end(); ++it) {
209  int stLayer = stationfromDet(it->second.front());
210  DTRecSegment4D bestDTSeg;
211  CSCSegment bestCSCSeg;
212  double bestDTEstimator = MAX_THR;
213  double bestCSCEstimator = MAX_THR;
214  vector<DTRecSegment4D> dtSegs = dtSegMap[it->first];
215  vector<CSCSegment> cscSegs = cscSegMap[it->first];
216 
217  // DT case: find the most compatible segment
218  TrajectoryStateOnSurface tsosDTlayer;
219  testDTstation(currentState, dtSegs, bestDTEstimator, bestDTSeg, tsosDTlayer);
220 
221  // CSC case: find the most compatible segment
222  TrajectoryStateOnSurface tsosCSClayer;
223  testCSCstation(currentState, cscSegs, bestCSCEstimator, bestCSCSeg, tsosCSClayer);
224 
225  // Decide whether to keep the layer or not
226  bool chosenLayer = chooseLayers(incompConLay, bestDTEstimator, bestDTSeg, tsosDTlayer, bestCSCEstimator, bestCSCSeg, tsosCSClayer);
227  fillDYTInfos(stLayer, chosenLayer, incompConLay, bestDTEstimator, bestCSCEstimator, bestDTSeg, bestCSCSeg);
228  }
229  //cout << "Number of used stations = " << nStationsUsed << endl;
230 }
231 
232 
233 //===> stationfromDet
235  if (det.subdetId() == MuonSubdetId::CSC) {
236  CSCDetId ch(det);
237  return ch.station();
238  }
239  if (det.subdetId() == MuonSubdetId::DT) {
240  DTChamberId ch(det);
241  return ch.station();
242  }
243  return 0;
244 }
245 
246 
247 //===> fillDYTInfos
248 void DynamicTruncation::fillDYTInfos(int const &st, bool const &chosenLayer, int &incompConLay,
249  double const &bestDTEstimator, double const &bestCSCEstimator,
250  DTRecSegment4D const &bestDTSeg, CSCSegment const &bestCSCSeg) {
251  if (chosenLayer) {
252  nStationsUsed++;
253  incompConLay = 0;
254  if (bestDTEstimator <= bestCSCEstimator) {
255  estimatorMap[st] = bestDTEstimator;
256  DetId id(bestDTSeg.chamberId());
257  idChamberMap[st] = id;
258  } else {
259  DetId id(bestCSCSeg.cscDetId());
260  idChamberMap[st] = id;
261  estimatorMap[st] = bestCSCEstimator;
262  }
263  usedStationMap[st] = true;
264  } else {
265  incompConLay++;
266  estimatorMap[st] = -1;
267  usedStationMap[st] = false;
268  }
269 }
270 
271 
272 //===> compatibleDets
273 void DynamicTruncation::compatibleDets(TrajectoryStateOnSurface &tsos, map<int, vector<DetId> > &detMap) {
274  MuonPatternRecoDumper dumper;
275  MeasurementEstimator *theEstimator = new Chi2MeasurementEstimator(1000, 1000);
276  vector<const DetLayer *> navLayers;
277  navLayers = navigation->compatibleLayers(*(currentState.freeState()), alongMomentum);
278  unsigned int ilayerCorrected = 0;
279  for ( unsigned int ilayer=0; ilayer<navLayers.size(); ilayer++ ) {
280  // Skip RPC layers
281  if (navLayers[ilayer]->subDetector() != GeomDetEnumerators::DT &&
282  navLayers[ilayer]->subDetector() != GeomDetEnumerators::CSC) continue;
283  ilayerCorrected++;
284  vector<DetLayer::DetWithState> comps = navLayers[ilayer]->compatibleDets(currentState, *propagatorCompatibleDet, *theEstimator);
285  //cout << comps.size() << " compatible Dets with " << navLayers[ilayer]->subDetector() << " Layer " << ilayer << " "
286  //<< dumper.dumpLayer(navLayers[ilayer]);
287  if (!comps.empty()) {
288  for ( unsigned int icomp=0; icomp<comps.size(); icomp++ ) {
289  DetId id(comps[icomp].first->geographicalId().rawId());
290  detMap[ilayerCorrected].push_back(id);
291  }
292  }
293  }
294  if (theEstimator) delete theEstimator;
295 }
296 
297 
298 //===> fillSegmentMaps
299 void DynamicTruncation::fillSegmentMaps( map<int, vector<DetId> > &compatibleIds,
300  map<int, vector<DTRecSegment4D> > &dtSegMap,
301  map<int, vector<CSCSegment> > &cscSegMap) {
302  for (map<int, vector<DetId> >::iterator it=compatibleIds.begin(); it!=compatibleIds.end(); ++it) {
303  vector<DetId> ids = compatibleIds[it->first];
304  for (unsigned j = 0; j < ids.size(); j++) {
305  if (ids[j].subdetId() == MuonSubdetId::CSC) {
306  CSCDetId ch(ids[j]);
307  vector<CSCSegment> tmp = getSegs->getCSCSegmentsInChamber(ch);
308  for (unsigned int k = 0; k < tmp.size(); k++) cscSegMap[it->first].push_back(tmp[k]);
309  }
310  if (ids[j].subdetId() == MuonSubdetId::DT) {
311  DTChamberId ch(ids[j]);
312  vector<DTRecSegment4D> tmp = getSegs->getDTSegmentsInChamber(ch);
313  for (unsigned int k = 0; k < tmp.size(); k++) dtSegMap[it->first].push_back(tmp[k]);
314  }
315  }
316  }
317 }
318 
319 
320 //===> testDTstation
321 void DynamicTruncation::testDTstation(TrajectoryStateOnSurface &startingState, vector<DTRecSegment4D> const &segments,
322  double &bestEstimator, DTRecSegment4D &bestSeg, TrajectoryStateOnSurface &tsosdt) {
323  if (segments.empty()) return;
324  for (unsigned int iSeg = 0; iSeg < segments.size(); iSeg++) {
325  DTChamberId chamber(segments[iSeg].chamberId());
326  if (!propagator->propagate(startingState, theG->idToDet(chamber)->surface()).isValid()) continue;
327  tsosdt = propagator->propagate(startingState, theG->idToDet(chamber)->surface());
328  //if (!tsosdt.isValid()) continue;
329  LocalError apeLoc;
330  if (useAPE) apeLoc = ErrorFrameTransformer().transform(dtApeMap.find(chamber)->second, theG->idToDet(chamber)->surface());
331  StateSegmentMatcher estim(tsosdt, segments[iSeg], apeLoc);
332  double estimator = estim.value();
333  //cout << "estimator DT = " << estimator << endl;
334  if (estimator >= bestEstimator) continue;
335  bestEstimator = estimator;
336  bestSeg = segments[iSeg];
337  }
338 }
339 
340 
341 //===> testCSCstation
342 void DynamicTruncation::testCSCstation(TrajectoryStateOnSurface &startingState, vector<CSCSegment> const &segments,
343  double &bestEstimator, CSCSegment &bestSeg, TrajectoryStateOnSurface &tsoscsc) {
344  if (segments.empty()) return;
345  for (unsigned int iSeg = 0; iSeg < segments.size(); iSeg++) {
346  CSCDetId chamber(segments[iSeg].cscDetId());
347  if (!propagator->propagate(startingState, theG->idToDet(chamber)->surface()).isValid()) continue;
348  tsoscsc = propagator->propagate(startingState, theG->idToDet(chamber)->surface());
349  //if (!tsoscsc.isValid()) continue;
350  LocalError apeLoc;
351  if (useAPE) apeLoc = ErrorFrameTransformer().transform(cscApeMap.find(chamber)->second, theG->idToDet(chamber)->surface());
352  StateSegmentMatcher estim(tsoscsc, segments[iSeg], apeLoc);
353  double estimator = estim.value();
354  //cout << "estimator CSC = " << estimator << endl;
355  if (estimator >= bestEstimator) continue;
356  bestEstimator = estimator;
357  bestSeg = segments[iSeg];
358  }
359 }
360 
361 
362 //===> useSegment
364  result.push_back(theMuonRecHitBuilder->build(&bestDTSeg));
365  if (doUpdateOfKFStates) updateWithDThits(currentState, bestDTSeg);
366  else currentState = tsosDT;
367 }
368 
369 
370 //===> useSegment
371 void DynamicTruncation::useSegment(CSCSegment const &bestCSCSeg, TrajectoryStateOnSurface const &tsosCSC) {
372  result.push_back(theMuonRecHitBuilder->build(&bestCSCSeg));
373  if (doUpdateOfKFStates) updateWithCSChits(currentState, bestCSCSeg);
374  else currentState = tsosCSC;
375 }
376 
377 
378 //===> preliminaryFit
379 void DynamicTruncation::preliminaryFit(map<int, vector<DetId> > compatibleIds, map<int, vector<DTRecSegment4D> > dtSegMap,
380  map<int, vector<CSCSegment> > cscSegMap) {
381  for (map<int, vector<DetId> >::iterator it=compatibleIds.begin(); it!=compatibleIds.end(); ++it) {
382  DTRecSegment4D bestDTSeg;
383  CSCSegment bestCSCSeg;
384  double bestDTEstimator = MAX_THR;
385  double bestCSCEstimator = MAX_THR;
386  double initThr = MAX_THR;
387  vector<DTRecSegment4D> dtSegs = dtSegMap[it->first];
388  vector<CSCSegment> cscSegs = cscSegMap[it->first];
389 
390  // DT case: find the most compatible segment
391  TrajectoryStateOnSurface tsosDTlayer;
392  testDTstation(prelFitState, dtSegs, bestDTEstimator, bestDTSeg, tsosDTlayer);
393 
394  // CSC case: find the most compatible segment
395  TrajectoryStateOnSurface tsosCSClayer;
396  testCSCstation(prelFitState, cscSegs, bestCSCEstimator, bestCSCSeg, tsosCSClayer);
397 
398  // Decide whether to keep the layer or not
399  if (bestDTEstimator == MAX_THR && bestCSCEstimator == MAX_THR) continue;
400  if (bestDTEstimator <= bestCSCEstimator) {
401  getThresholdFromCFG(initThr, DetId(bestDTSeg.chamberId()));
402  if (bestDTEstimator >= initThr) continue;
403  prelFitMeas.push_back(theMuonRecHitBuilder->build(&bestDTSeg));
404  auto aSegRH = prelFitMeas.back();
405  auto uRes = updatorHandle->update(tsosDTlayer, *aSegRH);
406  if (uRes.isValid()){
407  prelFitState = uRes;
408  } else {
409  prelFitMeas.pop_back();
410  }
411  } else {
412  getThresholdFromCFG(initThr, DetId(bestCSCSeg.cscDetId()));
413  if (bestCSCEstimator >= initThr) continue;
414  prelFitMeas.push_back(theMuonRecHitBuilder->build(&bestCSCSeg));
415  auto aSegRH = prelFitMeas.back();
416  auto uRes = updatorHandle->update(tsosCSClayer, *aSegRH);
417  if (uRes.isValid()){
418  prelFitState = uRes;
419  } else {
420  prelFitMeas.pop_back();
421  }
422  }
423  }
424  if (!prelFitMeas.empty()) prelFitMeas.pop_back();
425  for (auto imrh = prelFitMeas.rbegin(); imrh != prelFitMeas.rend(); ++imrh) {
426  DetId id = (*imrh)->geographicalId();
427  TrajectoryStateOnSurface tmp = propagatorPF->propagate(prelFitState, theG->idToDet(id)->surface());
428  if (tmp.isValid()) prelFitState = tmp;
429  }
430  muonPTest = prelFitState.globalMomentum().perp();
431  muonETAest = prelFitState.globalMomentum().eta();
432 }
433 
434 
435 //===> chooseLayers
436 bool DynamicTruncation::chooseLayers(int &incompLayers, double const &bestDTEstimator, DTRecSegment4D const &bestDTSeg, TrajectoryStateOnSurface const &tsosDT,
437  double const &bestCSCEstimator, CSCSegment const &bestCSCSeg, TrajectoryStateOnSurface const &tsosCSC) {
438  double initThr = MAX_THR;
439  if (bestDTEstimator == MAX_THR && bestCSCEstimator == MAX_THR) return false;
440  if (bestDTEstimator <= bestCSCEstimator) {
441  // Get threshold for the chamber
442  if (useDBforThr) getThresholdFromDB(initThr, DetId(bestDTSeg.chamberId()));
443  else getThresholdFromCFG(initThr, DetId(bestDTSeg.chamberId()));
444  if (DYTselector == 0 || (DYTselector == 1 && bestDTEstimator < initThr) ||
445  (DYTselector == 2 && incompLayers < 2 && bestDTEstimator < initThr)) {
446  useSegment(bestDTSeg, tsosDT);
447  return true;
448  }
449  } else {
450  // Get threshold for the chamber
451  if (useDBforThr) getThresholdFromDB(initThr, DetId(bestCSCSeg.cscDetId()));
452  else getThresholdFromCFG(initThr, DetId(bestCSCSeg.cscDetId()));
453  if (DYTselector == 0 || (DYTselector == 1 && bestCSCEstimator < initThr) ||
454  (DYTselector == 2 && incompLayers < 2 && bestCSCEstimator < initThr)) {
455  useSegment(bestCSCSeg, tsosCSC);
456  return true;
457  }
458  }
459  return false;
460 }
461 
462 
463 //===> getThresholdFromDB
464 void DynamicTruncation::getThresholdFromDB(double& thr, DetId const& id) {
465  vector<DYTThrObject::DytThrStruct> thrvector = dytThresholds->thrsVec;
466  for (vector<DYTThrObject::DytThrStruct>::const_iterator it = thrvector.begin(); it != thrvector.end(); it++) {
468  if (obj.id == id) {
469  thr = obj.thr;
470  break;
471  }
472  }
473  if (useParametrizedThr) correctThrByPAndEta(thr);
474 }
475 
476 
477 //===> correctThrByPAndEta
479 
480  auto parametricThreshold = [this]{
481  double thr50 = this->parameters[this->region].at(0);
482  double p0 = this->parameters[this->region].at(1);
483  double p1 = this->parameters[this->region].at(2);
484  return thr50 * ( 1 + p0*p_reco + std::pow( this->p_reco, p1));
485  };
486 
487  thr = parametricThreshold();
488 }
489 
491 
492  float absEta = std::abs(eta_reco);
493  // Defaul value for muons with abs(eta) > 2.4
495 
496  if( absEta <= 0.8 ) region = dyt_utils::etaRegion::eta0p8;
497  else if( absEta <= 1.2 ) region = dyt_utils::etaRegion::eta1p2;
498  else if( absEta <= 2.0 ) region = dyt_utils::etaRegion::eta2p0;
499  else if( absEta <= 2.2 ) region = dyt_utils::etaRegion::eta2p2;
500  else if( absEta <= 2.4 ) region = dyt_utils::etaRegion::eta2p4;
501 
502 }
503 
504 
505 //===> getThresholdFromCFG
506 void DynamicTruncation::getThresholdFromCFG(double& thr, DetId const& id) {
507  if (id.subdetId() == MuonSubdetId::DT) {
508  thr = Thrs[0];
509  }
510  if (id.subdetId() == MuonSubdetId::CSC) {
511  thr = Thrs[1];
512  }
513 
514  if (useParametrizedThr) correctThrByPAndEta(thr);
515 }
516 
517 
518 //===> sort
520  unsigned int i=0;
521  unsigned int j=0;
522  ConstRecHitContainer::size_type n = recHits.size();
523  for(i=1; i<n; ++i)
524  for(j=n-1; j>=i; --j)
525  if(recHits[j-1]->globalPosition().mag() > recHits[j]->globalPosition().mag()) swap (recHits[j-1],recHits[j]);
526 }
T getParameter(std::string const &) const
static GlobalError transform(const LocalError &le, const Surface &surf)
void compatibleDets(TrajectoryStateOnSurface &, std::map< int, std::vector< DetId > > &)
T perp() const
Definition: PV3DBase.h:72
CSCDetId cscDetId() const
Definition: CSCSegment.h:69
ConstRecHitPointer const & recHit() const
void sort(ConstRecHitContainer &)
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
std::vector< const TrackingRecHit * > recHits() const override
Access to component RecHits (if any)
void fillSegmentMaps(std::map< int, std::vector< DetId > > &, std::map< int, std::vector< DTRecSegment4D > > &, std::map< int, std::vector< CSCSegment > > &)
void useSegment(DTRecSegment4D const &, TrajectoryStateOnSurface const &)
GlobalPoint globalPosition() const
void getThresholdFromDB(double &, DetId const &)
virtual DTChamberId chamberId() const
The (specific) DetId of the chamber on which the segment resides.
void updateWithCSChits(TrajectoryStateOnSurface &, CSCSegment const &)
void testDTstation(TrajectoryStateOnSurface &, std::vector< DTRecSegment4D > const &, double &, DTRecSegment4D &, TrajectoryStateOnSurface &)
uint16_t size_type
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:42
void swap(Association< C > &lhs, Association< C > &rhs)
Definition: Association.h:116
DataContainer const & measurements() const
Definition: Trajectory.h:196
static const int CSC
Definition: MuonSubdetId.h:13
const SurfaceType & surface() const
T mag() const
Definition: PV3DBase.h:67
TransientTrackingRecHit::ConstRecHitContainer ConstRecHitContainer
void update(const LocalTrajectoryParameters &p, const SurfaceType &aSurface, const MagneticField *field, SurfaceSide side=SurfaceSideDefinition::atCenterOfSurface)
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:41
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void correctThrByPAndEta(double &)
void testCSCstation(TrajectoryStateOnSurface &, std::vector< CSCSegment > const &, double &, CSCSegment &, TrajectoryStateOnSurface &)
const bool isValidThdDB()
Definition: ThrParameters.h:28
TransientTrackingRecHit::ConstRecHitContainer filter(const Trajectory &)
const std::vector< CSCRecHit2D > & specificRecHits() const
Definition: CSCSegment.h:65
#define MAX_THR
TrajectoryStateOnSurface const & forwardPredictedState() const
Access to forward predicted state (from fitter or builder)
std::vector< ConstRecHitPointer > ConstRecHitContainer
int k[5][pyjets_maxn]
void preliminaryFit(std::map< int, std::vector< DetId > >, std::map< int, std::vector< DTRecSegment4D > >, std::map< int, std::vector< CSCSegment > >)
void setThrsMap(const edm::ParameterSet &)
Definition: DetId.h:18
DynamicTruncation(const edm::Event &, const MuonServiceProxy &)
void update(TrajectoryStateOnSurface &, ConstRecHitPointer)
void fillDYTInfos(int const &, bool const &, int &, double const &, double const &, DTRecSegment4D const &, CSCSegment const &)
TransientTrackingRecHit::ConstRecHitPointer ConstRecHitPointer
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
fixed size matrix
HLT enums.
double p1[4]
Definition: TauolaWrapper.h:89
GlobalVector globalMomentum() const
#define update(a, b)
int station() const
Definition: CSCDetId.h:86
static const int DT
Definition: MuonSubdetId.h:12
bool chooseLayers(int &, double const &, DTRecSegment4D const &, TrajectoryStateOnSurface const &, double const &, CSCSegment const &, TrajectoryStateOnSurface const &)
int stationfromDet(DetId const &)
DetId geographicalId() const
static const std::map< etaRegion, std::string > etaRegionStr
int station() const
Return the station number.
Definition: DTChamberId.h:51
T first(std::pair< T, U > const &p)
void updateWithDThits(TrajectoryStateOnSurface &, DTRecSegment4D const &)
void setThr(const std::vector< int > &)
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
Definition: event.py:1
void getThresholdFromCFG(double &, DetId const &)
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:39