CMS 3D CMS Logo

ConvBremSeedProducer.cc
Go to the documentation of this file.
49 
50 #include "TMath.h"
51 
52 #include <memory>
53 
55 public:
56  explicit ConvBremSeedProducer(const edm::ParameterSet&);
57 
58  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
59 
60 private:
61  void beginRun(const edm::Run&, const edm::EventSetup&) override;
62  void produce(edm::Event&, const edm::EventSetup&) override;
63  void endRun(const edm::Run&, const edm::EventSetup&) override;
64  void initializeLayerMap();
65  std::vector<const DetLayer*> theLayerMap;
67  const ParticlePropagator& pp,
68  const MagneticField* field) const;
69  const DetLayer* detLayer(const TrackerLayer& layer, float zpos) const;
70 
71  bool isGsfTrack(const reco::Track&, const TrackingRecHit*);
72 
73  int GoodCluster(const BaseParticlePropagator& bpg,
74  const reco::PFClusterCollection& pfc,
75  float minep,
76  bool sec = false);
77 
78  std::vector<bool> sharedHits(const std::vector<std::pair<TrajectorySeed, std::pair<GlobalVector, float> > >&);
79 
89  std::vector<const DetLayer*> layerMap_;
93 
94  // Event setup tokens
102 };
103 
106 
108  // convBremSeeds
110  desc.add<edm::InputTag>("pixelRecHits", edm::InputTag("gsPixelRecHits"));
111  desc.add<edm::InputTag>("matchedrecHits", edm::InputTag("gsStripRecHits", "matchedRecHit"));
112  desc.add<std::string>("TTRHBuilder", "WithTrackAngle");
113  desc.add<edm::InputTag>("rphirecHits", edm::InputTag("gsStripRecHits", "rphiRecHit"));
114  desc.add<edm::InputTag>("PFClusters", edm::InputTag("particleFlowClusterECAL"));
115  desc.add<edm::InputTag>("PFRecTrackLabel", edm::InputTag("pfTrackElec"));
116 }
117 
118 using namespace edm;
119 using namespace std;
120 using namespace reco;
121 
123  : conf_(iConfig),
124  fieldMap_(nullptr),
125  layerMap_(56, static_cast<const DetLayer*>(nullptr)),
126  negLayerOffset_(27),
127  magFieldToken_(esConsumes()),
128  geomSearchTrackerToken_(esConsumes<edm::Transition::BeginRun>()),
129  geometryToken_(esConsumes<edm::Transition::BeginRun>()),
130  trackerToken_(esConsumes<edm::Transition::BeginRun>()),
131  magFieldToken_beginRun_(esConsumes<edm::Transition::BeginRun>()),
132  magFieldMapToken_(esConsumes<edm::Transition::BeginRun>()),
133  hitBuilderToken_(
134  esConsumes<edm::Transition::BeginRun>(edm::ESInputTag("", conf_.getParameter<string>("TTRHBuilder")))) {
135  produces<ConvBremSeedCollection>();
136 }
137 
139  LogDebug("ConvBremSeedProducerProducer") << "START event: " << iEvent.id().event() << " in run " << iEvent.id().run();
140 
141  constexpr float pfmass = 0.0005;
142 
144 
147  iEvent.getByLabel(conf_.getParameter<InputTag>("PFClusters"), PfC);
148  const PFClusterCollection& PPP = *(PfC.product());
149 
152  iEvent.getByLabel(conf_.getParameter<InputTag>("pixelRecHits"), pixelHits);
153 
156  iEvent.getByLabel(conf_.getParameter<InputTag>("rphirecHits"), rphirecHits);
158  iEvent.getByLabel(conf_.getParameter<InputTag>("matchedrecHits"), matchedrecHits);
159 
160  //GSFPFRECTRACKS
161  Handle<GsfPFRecTrackCollection> thePfRecTrackCollection;
162  iEvent.getByLabel(conf_.getParameter<InputTag>("PFRecTrackLabel"), thePfRecTrackCollection);
163  const GsfPFRecTrackCollection PfRTkColl = *(thePfRecTrackCollection.product());
164 
166  auto output = std::make_unique<ConvBremSeedCollection>();
167 
169  vector<pair<TrajectorySeed, pair<GlobalVector, float> > > unclean;
170  //TRIPLET OF MODULES TO BE USED FOR SEEDING
171  vector<vector<long int> > tripl;
172  //LAYER MAP
174 
176 
177  for (unsigned int ipft = 0; ipft < PfRTkColl.size(); ipft++) {
178  GsfPFRecTrackRef pft(thePfRecTrackCollection, ipft);
179  LogDebug("ConvBremSeedProducerProducer") << "NEW GsfPFRecTRACK ";
180  float eta_br = 0;
181  unclean.clear();
182  tripl.clear();
183  vector<int> gc;
184  auto const& gsfRecHits = *pft->gsfTrackRef();
185  float pfoutenergy = sqrt((pfmass * pfmass) + pft->gsfTrackRef()->outerMomentum().Mag2());
186  XYZTLorentzVector mom = XYZTLorentzVector(pft->gsfTrackRef()->outerMomentum().x(),
187  pft->gsfTrackRef()->outerMomentum().y(),
188  pft->gsfTrackRef()->outerMomentum().z(),
189  pfoutenergy);
190  XYZTLorentzVector pos = XYZTLorentzVector(pft->gsfTrackRef()->outerPosition().x(),
191  pft->gsfTrackRef()->outerPosition().y(),
192  pft->gsfTrackRef()->outerPosition().z(),
193  0.);
194  BaseParticlePropagator theOutParticle =
195  BaseParticlePropagator(RawParticle(mom, pos, pft->gsfTrackRef()->charge()), 0, 0, B_.z());
196 
198  gc.push_back(GoodCluster(theOutParticle, PPP, 0.5));
199 
200  vector<PFBrem> brem = (*pft).PFRecBrem();
201  vector<PFBrem>::iterator ib = brem.begin();
202  vector<PFBrem>::iterator ib_end = brem.end();
203  LogDebug("ConvBremSeedProducerProducer") << "NUMBER OF BREMS " << brem.size();
204 
206  for (; ib != ib_end; ++ib) {
207  XYZTLorentzVector mom = pft->trajectoryPoint(ib->indTrajPoint()).momentum();
208  XYZTLorentzVector pos = XYZTLorentzVector(pft->trajectoryPoint(ib->indTrajPoint()).position().x(),
209  pft->trajectoryPoint(ib->indTrajPoint()).position().y(),
210  pft->trajectoryPoint(ib->indTrajPoint()).position().z(),
211  0);
212 
214  if (pos.Rho() > 80)
215  continue;
216  if ((pos.Rho() > 5) && (fabs(ib->SigmaDeltaP() / ib->DeltaP()) > 3))
217  continue;
218  if (fabs(ib->DeltaP()) < 3)
219  continue;
220  eta_br = mom.eta();
221  vector<vector<long int> > Idd;
222 
223  BaseParticlePropagator p(RawParticle(mom, pos, 0.), 0, 0, B_.z());
224  gc.push_back(GoodCluster(p, PPP, 0.2));
225 
226  ParticlePropagator PP(p, fieldMap_, nullptr);
227 
229  list<TrackerLayer>::const_iterator cyliter = geometry_->cylinderBegin();
230  for (; cyliter != geometry_->cylinderEnd(); ++cyliter) {
232  if (!(cyliter->sensitive()))
233  continue;
234  PP.setPropagationConditions(*cyliter);
235  PP.propagate();
236  if (PP.getSuccess() == 0)
237  continue;
238 
241  AnalyticalPropagator alongProp(&mf, anyDirection);
243  const DetLayer* tkLayer = detLayer(*cyliter, PP.particle().Z());
244  if (&(*tkLayer) == nullptr)
245  continue;
246  TrajectoryStateOnSurface trajState = makeTrajectoryState(tkLayer, PP, &mf);
247 
248  auto compat = tkLayer->compatibleDets(trajState, alongProp, est);
249  vector<long int> temp;
250  if (compat.empty())
251  continue;
252 
253  for (auto i = compat.begin(); i != compat.end(); i++) {
254  long int detid = i->first->geographicalId().rawId();
255 
257  auto DetMatch = (rphirecHits.product())->find((detid));
258  auto MDetMatch = (matchedrecHits.product())->find((detid));
259 
260  long int DetID = (DetMatch != rphirecHits->end()) ? detid : 0;
261 
262  if ((MDetMatch != matchedrecHits->end()) && !MDetMatch->empty()) {
263  long int pii = MDetMatch->begin()->monoId();
264  auto CDetMatch = (rphirecHits.product())->find((pii));
265  DetID = (CDetMatch != rphirecHits->end()) ? pii : 0;
266  }
267 
268  temp.push_back(DetID);
269 
270  } else {
271  auto DetMatch = (pixelHits.product())->find((detid));
272  long int DetID = (DetMatch != pixelHits->end()) ? detid : 0;
273  temp.push_back(DetID);
274  }
275  }
276 
277  Idd.push_back(temp);
278 
279  } //END TRACKER LAYER LOOP
280  if (Idd.size() < 2)
281  continue;
282 
284  for (unsigned int i = 0; i < Idd.size() - 2; i++) {
285  for (unsigned int i1 = 0; i1 < Idd[i].size(); i1++) {
286  for (unsigned int i2 = 0; i2 < Idd[i + 1].size(); i2++) {
287  for (unsigned int i3 = 0; i3 < Idd[i + 2].size(); i3++) {
288  if ((Idd[i][i1] != 0) && (Idd[i + 1][i2] != 0) && (Idd[i + 2][i3] != 0)) {
289  vector<long int> tmp;
290  tmp.push_back(Idd[i][i1]);
291  tmp.push_back(Idd[i + 1][i2]);
292  tmp.push_back(Idd[i + 2][i3]);
293 
294  bool newTrip = true;
295  for (unsigned int iv = 0; iv < tripl.size(); iv++) {
296  if ((tripl[iv][0] == tmp[0]) && (tripl[iv][1] == tmp[1]) && (tripl[iv][2] == tmp[2]))
297  newTrip = false;
298  }
299  if (newTrip) {
300  tripl.push_back(tmp);
301  }
302  }
303  }
304  }
305  }
306  }
307  } //END BREM LOOP
308 
309  float sineta_brem = sinh(eta_br);
310 
311  //OUTPUT COLLECTION
312  auto bfield = iSetup.getHandle(magFieldToken_);
313  float nomField = bfield->nominalValue();
314 
316  OwnVector<TrackingRecHit> loc_hits;
317  for (unsigned int i = 0; i < tripl.size(); i++) {
318  auto DetMatch1 = (rphirecHits.product())->find(tripl[i][0]);
319  auto DetMatch2 = (rphirecHits.product())->find(tripl[i][1]);
320  auto DetMatch3 = (rphirecHits.product())->find(tripl[i][2]);
321  if ((DetMatch1 == rphirecHits->end()) || (DetMatch2 == rphirecHits->end()) || (DetMatch3 == rphirecHits->end()))
322  continue;
323  auto DetSet1 = *DetMatch1;
324  auto DetSet2 = *DetMatch2;
325  auto DetSet3 = *DetMatch3;
326 
327  for (auto it1 = DetSet1.begin(); it1 != DetSet1.end(); ++it1) {
328  GlobalPoint gp1 = tracker_->idToDet(tripl[i][0])->surface().toGlobal(it1->localPosition());
329 
330  bool tak1 = isGsfTrack(gsfRecHits, &(*it1));
331 
332  for (auto it2 = DetSet2.begin(); it2 != DetSet2.end(); ++it2) {
333  GlobalPoint gp2 = tracker_->idToDet(tripl[i][1])->surface().toGlobal(it2->localPosition());
334  bool tak2 = isGsfTrack(gsfRecHits, &(*it2));
335 
336  for (auto it3 = DetSet3.begin(); it3 != DetSet3.end(); ++it3) {
337  // ips++;
338  GlobalPoint gp3 = tracker_->idToDet(tripl[i][2])->surface().toGlobal(it3->localPosition());
339  bool tak3 = isGsfTrack(gsfRecHits, &(*it3));
340 
341  FastHelix helix(gp3, gp2, gp1, nomField, &*bfield);
342  GlobalVector gv = helix.stateAtVertex().momentum();
343  GlobalVector gv_corr(gv.x(), gv.y(), gv.perp() * sineta_brem);
344  float ene = sqrt(gv_corr.mag2() + (pfmass * pfmass));
345 
346  GlobalPoint gp = helix.stateAtVertex().position();
347  float ch = helix.stateAtVertex().charge();
348 
349  XYZTLorentzVector mom = XYZTLorentzVector(gv.x(), gv.y(), gv_corr.z(), ene);
350  XYZTLorentzVector pos = XYZTLorentzVector(gp.x(), gp.y(), gp.z(), 0.);
351  BaseParticlePropagator theOutParticle(RawParticle(mom, pos, ch), 0, 0, B_.z());
352  int bgc = GoodCluster(theOutParticle, PPP, 0.3, true);
353 
354  if (gv.perp() < 0.5)
355  continue;
356 
357  if (tak1 + tak2 + tak3 > 2)
358  continue;
359 
360  if (bgc == -1)
361  continue;
362  bool clTak = false;
363  for (unsigned int igcc = 0; igcc < gc.size(); igcc++) {
364  if (clTak)
365  continue;
366  if (bgc == gc[igcc])
367  clTak = true;
368  }
369  if (clTak)
370  continue;
371 
372  GlobalTrajectoryParameters Gtp(gp1, gv, int(ch), &(*magfield_));
373  glob_hits.clear();
374  loc_hits.clear();
375  glob_hits.push_back(hitBuilder_->build(it1->clone()));
376  glob_hits.push_back(hitBuilder_->build(it2->clone()));
377  glob_hits.push_back(hitBuilder_->build(it3->clone()));
378 
380 
382  TrajectoryStateOnSurface updatedState;
383 
384  for (int ih = 0; ih < 3; ih++) {
386  (ih == 0) ? propagator_->propagate(CSeed, tracker_->idToDet(tripl[i][ih])->surface())
387  : propagator_->propagate(updatedState, tracker_->idToDet(tripl[i][ih])->surface());
388 
389  if (!state.isValid()) {
390  ih = 3;
391  continue;
392  }
393 
394  updatedState = kfUpdator_->update(state, *glob_hits[ih]);
395  loc_hits.push_back(glob_hits[ih]->hit()->clone());
396  if (ih == 2) {
397  PTrajectoryStateOnDet const& PTraj =
398  trajectoryStateTransform::persistentState(updatedState, tripl[i][2]);
399  // output->push_back(Trajectoryseed(PTraj,loc_hits,alongMomentum));
400  unclean.push_back(make_pair(TrajectorySeed(PTraj, loc_hits, alongMomentum), make_pair(gv_corr, ch)));
401  }
402  // }
403  }
404  }
405  }
406  }
407  }
408  vector<bool> inPhot = sharedHits(unclean);
409  for (unsigned int iu = 0; iu < unclean.size(); iu++) {
410  if (inPhot[iu])
411  output->push_back(ConvBremSeed(unclean[iu].first, pft));
412  }
413 
414  } //END GSF TRACK COLLECTION LOOP
415  LogDebug("ConvBremSeedProducerProducer") << "END";
416  iEvent.put(std::move(output));
417 }
418 
421 
422  geometry_ = &iSetup.getData(geometryToken_);
423 
424  tracker_ = &iSetup.getData(trackerToken_);
425 
427  B_ = magfield_->inTesla(GlobalPoint(0, 0, 0));
428 
430 
432 
434  kfUpdator_ = new KFUpdator();
435 }
436 
438  delete propagator_;
439  delete kfUpdator_;
440 }
441 
443  // These are the BoundSurface&, the BoundDisk* and the BoundCylinder* for that layer
444  // const BoundSurface& theSurface = layer.surface();
445  // BoundDisk* theDisk = layer.disk(); // non zero for endcaps
446  // BoundCylinder* theCylinder = layer.cylinder(); // non zero for barrel
447  // int theLayer = layer.layerNumber(); // 1->3 PixB, 4->5 PixD,
448  // // 6->9 TIB, 10->12 TID,
449  // // 13->18 TOB, 19->27 TEC
450 
453 
454  const std::vector<const BarrelDetLayer*>& barrelLayers = geomSearchTracker_->barrelLayers();
455  LogDebug("FastTracker") << "Barrel DetLayer dump: ";
456  for (auto bl = barrelLayers.begin(); bl != barrelLayers.end(); ++bl) {
457  LogDebug("FastTracker") << "radius " << (**bl).specificSurface().radius();
458  }
459 
460  const std::vector<const ForwardDetLayer*>& posForwardLayers = geomSearchTracker_->posForwardLayers();
461  LogDebug("FastTracker") << "Positive Forward DetLayer dump: ";
462  for (auto fl = posForwardLayers.begin(); fl != posForwardLayers.end(); ++fl) {
463  LogDebug("FastTracker") << "Z pos " << (**fl).surface().position().z() << " radii "
464  << (**fl).specificSurface().innerRadius() << ", " << (**fl).specificSurface().outerRadius();
465  }
466 
467  const float rTolerance = 1.5;
468  const float zTolerance = 3.;
469 
470  LogDebug("FastTracker") << "Dump of TrackerInteractionGeometry cylinders:";
471  for (std::list<TrackerLayer>::const_iterator i = geometry_->cylinderBegin(); i != geometry_->cylinderEnd(); ++i) {
472  const BoundCylinder* cyl = i->cylinder();
473  const BoundDisk* disk = i->disk();
474 
475  LogDebug("FastTracker") << "Famos Layer no " << i->layerNumber() << " is sensitive? " << i->sensitive() << " pos "
476  << i->surface().position();
477  if (!i->sensitive())
478  continue;
479 
480  if (cyl != nullptr) {
481  LogDebug("FastTracker") << " cylinder radius " << cyl->radius();
482  bool found = false;
483 
484  for (auto bl = barrelLayers.begin(); bl != barrelLayers.end(); ++bl) {
485  if (fabs(cyl->radius() - (**bl).specificSurface().radius()) < rTolerance) {
486  layerMap_[i->layerNumber()] = *bl;
487  found = true;
488  LogDebug("FastTracker") << "Corresponding DetLayer found with radius " << (**bl).specificSurface().radius();
489 
490  break;
491  }
492  }
493  if (!found) {
494  LogError("FastTracker") << "FAILED to find a corresponding DetLayer!";
495  }
496  } else {
497  LogDebug("FastTracker") << " disk radii " << disk->innerRadius() << ", " << disk->outerRadius();
498 
499  bool found = false;
500 
501  for (auto fl = posForwardLayers.begin(); fl != posForwardLayers.end(); ++fl) {
502  if (fabs(disk->position().z() - (**fl).surface().position().z()) < zTolerance) {
503  layerMap_[i->layerNumber()] = *fl;
504  found = true;
505  LogDebug("FastTracker") << "Corresponding DetLayer found with Z pos " << (**fl).surface().position().z()
506  << " and radii " << (**fl).specificSurface().innerRadius() << ", "
507  << (**fl).specificSurface().outerRadius();
508  break;
509  }
510  }
511  if (!found) {
512  LogError("FastTracker") << "FAILED to find a corresponding DetLayer!";
513  }
514  }
515  }
516 }
517 const DetLayer* ConvBremSeedProducer::detLayer(const TrackerLayer& layer, float zpos) const {
518  if (zpos > 0 || !layer.forward())
519  return layerMap_[layer.layerNumber()];
520  else
521  return layerMap_[layer.layerNumber() + negLayerOffset_];
522 }
523 
525  const ParticlePropagator& pp,
526  const MagneticField* field) const {
527  GlobalPoint pos(pp.particle().X(), pp.particle().Y(), pp.particle().Z());
528  GlobalVector mom(pp.particle().Px(), pp.particle().Py(), pp.particle().Pz());
529 
530  auto plane = layer->surface().tangentPlane(pos);
531 
532  return TrajectoryStateOnSurface(GlobalTrajectoryParameters(pos, mom, TrackCharge(pp.particle().charge()), field),
533  *plane);
534 }
536  bool istaken = false;
537  for (auto const& hit : tkv.recHits()) {
538  if (istaken || !hit->isValid())
539  continue;
540  istaken = hit->sharesInput(h, TrackingRecHit::all);
541  }
542  return istaken;
543 }
544 vector<bool> ConvBremSeedProducer::sharedHits(const vector<pair<TrajectorySeed, pair<GlobalVector, float> > >& unclean) {
545  vector<bool> goodseed;
546  goodseed.clear();
547  if (unclean.size() < 2) {
548  for (unsigned int i = 0; i < unclean.size(); i++)
549  goodseed.push_back(true);
550  } else {
551  for (unsigned int i = 0; i < unclean.size(); i++)
552  goodseed.push_back(true);
553 
554  for (unsigned int iu = 0; iu < unclean.size() - 1; iu++) {
555  if (!goodseed[iu])
556  continue;
557  for (unsigned int iu2 = iu + 1; iu2 < unclean.size(); iu2++) {
558  if (!goodseed[iu])
559  continue;
560  if (!goodseed[iu2])
561  continue;
562  // if (unclean[iu].second.second *unclean[iu2].second.second >0)continue;
563 
564  unsigned int shar = 0;
565  for (auto const& sh : unclean[iu].first.recHits()) {
566  for (auto const& sh2 : unclean[iu2].first.recHits()) {
567  if (sh.sharesInput(&sh2, TrackingRecHit::all))
568 
569  shar++;
570  }
571  }
572  if (shar >= 2) {
573  if (unclean[iu].second.first.perp() < unclean[iu2].second.first.perp())
574  goodseed[iu] = false;
575  else
576  goodseed[iu2] = false;
577  }
578  }
579  }
580  }
581  return goodseed;
582 }
583 
585  const PFClusterCollection& pfc,
586  float minep,
587  bool sec) {
588  BaseParticlePropagator bpg = ubpg;
589  bpg.propagateToEcalEntrance(false);
590  float dr = 1000;
591  float de = 1000;
592  float df = 1000;
593  int ibest = -1;
594 
595  if (bpg.getSuccess() != 0) {
596  for (unsigned int i = 0; i < pfc.size(); i++) {
597  float tmp_ep = pfc[i].energy() / bpg.particle().momentum().e();
598  float tmp_phi = fabs(pfc[i].position().phi() - bpg.particle().vertex().phi());
599  if (tmp_phi > TMath::TwoPi())
600  tmp_phi -= TMath::TwoPi();
601  float tmp_eta = fabs(pfc[i].position().eta() - bpg.particle().vertex().eta());
602  float tmp_dr = sqrt(pow(tmp_phi, 2) + pow(tmp_eta, 2));
603  bool isBet = (tmp_dr < dr);
604  if (sec)
605  isBet = (tmp_phi < df);
606  if ((isBet) && (tmp_ep > minep) && (tmp_ep < 10)) {
607  dr = tmp_dr;
608  de = tmp_eta;
609  df = tmp_phi;
610  ibest = i;
611  }
612  }
613  bool isBad = (dr > 0.1);
614  if (sec)
615  isBad = ((df > 0.25) || (de > 0.5));
616 
617  if (isBad)
618  ibest = -1;
619  }
620  return ibest;
621 }
const double TwoPi
const MagneticField * magfield_
const GeometricSearchTracker * geomSearchTracker_
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
const TransientTrackingRecHitBuilder * hitBuilder_
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
ConvBremSeedProducer(const edm::ParameterSet &)
const edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > trackerToken_
T perp() const
Definition: PV3DBase.h:69
GlobalTrajectoryParameters stateAtVertex() const
Definition: FastHelix.h:59
virtual std::vector< DetWithState > compatibleDets(const TrajectoryStateOnSurface &startingState, const Propagator &prop, const MeasurementEstimator &est) const
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magFieldToken_
TrajectoryStateOnSurface makeTrajectoryState(const DetLayer *layer, const ParticlePropagator &pp, const MagneticField *field) const
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
const edm::ESGetToken< TransientTrackingRecHitBuilder, TransientRecHitRecord > hitBuilderToken_
virtual GlobalVector inTesla(const GlobalPoint &gp) const =0
Field value ad specified global point, in Tesla.
double Z() const
z of vertex
Definition: RawParticle.h:288
ROOT::Math::SMatrixIdentity AlgebraicMatrixID
const edm::ESGetToken< TrackerInteractionGeometry, TrackerInteractionGeometryRecord > geometryToken_
T const * product() const
Definition: Handle.h:70
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
PTrajectoryStateOnDet persistentState(const TrajectoryStateOnSurface &ts, unsigned int detid)
const PropagatorWithMaterial * propagator_
virtual SubDetector subDetector() const =0
The type of detector (PixelBarrel, PixelEndcap, TIB, TOB, TID, TEC, CSC, DT, RPCBarrel, RPCEndcap)
const XYZTLorentzVector & momentum() const
the momentum fourvector
Definition: RawParticle.h:321
Log< level::Error, false > LogError
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
const KFUpdator * kfUpdator_
const TrackerInteractionGeometry * geometry_
const MagneticFieldMap * fieldMap_
U second(std::pair< T, U > const &p)
bool isGsfTrack(const reco::Track &, const TrackingRecHit *)
auto recHits() const
Access to reconstructed hits on the track.
Definition: Track.h:85
int TrackCharge
Definition: TrackCharge.h:4
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
void push_back(D *&d)
Definition: OwnVector.h:326
int iEvent
Definition: GenABIO.cc:224
math::XYZVector B_
B field.
TrajectoryStateOnSurface update(const TrajectoryStateOnSurface &, const TrackingRecHit &) const override
Definition: KFUpdator.cc:177
std::vector< GsfPFRecTrack > GsfPFRecTrackCollection
collection of GsfPFRecTrack objects
void produce(edm::Event &, const edm::EventSetup &) override
T sqrt(T t)
Definition: SSEVec.h:19
virtual RecHitPointer build(const TrackingRecHit *p) const =0
build a tracking rechit from an existing rechit
void clear()
Definition: OwnVector.h:481
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void setPropagationConditions(const TrackerLayer &, bool firstLoop=true)
Transition
Definition: Transition.h:12
RawParticle const & particle() const
The particle being propagated.
std::list< TrackerLayer >::const_iterator cylinderEnd() const
Returns the last pointer in the cylinder list.
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::vector< ForwardDetLayer const * > const & posForwardLayers() const
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magFieldToken_beginRun_
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
const TrackerGeomDet * idToDet(DetId) const override
bool propagateToEcalEntrance(bool first=true)
std::vector< const DetLayer * > layerMap_
std::vector< ConstRecHitPointer > ConstRecHitContainer
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:79
std::vector< bool > sharedHits(const std::vector< std::pair< TrajectorySeed, std::pair< GlobalVector, float > > > &)
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
std::list< TrackerLayer >::const_iterator cylinderBegin() const
Returns the first pointer in the cylinder list.
int getSuccess() const
Has propagation been performed and was barrel or endcap reached ?
const edm::ESGetToken< MagneticFieldMap, MagneticFieldMapRecord > magFieldMapToken_
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
std::vector< BarrelDetLayer const * > const & barrelLayers() const
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
TEveGeoShape * clone(const TEveElement *element, TEveElement *parent)
Definition: eve_macros.cc:135
void beginRun(const edm::Run &, const edm::EventSetup &) override
TrajectoryStateOnSurface propagate(STA const &state, SUR const &surface) const
Definition: Propagator.h:50
void endRun(const edm::Run &, const edm::EventSetup &) override
const DetLayer * detLayer(const TrackerLayer &layer, float zpos) const
fixed size matrix
HLT enums.
std::vector< const DetLayer * > theLayerMap
static int position[264][3]
Definition: ReadPGInfo.cc:289
std::vector< PFCluster > PFClusterCollection
collection of PFCluster objects
Definition: PFClusterFwd.h:9
Definition: output.py:1
double getMagneticField() const
Get the magnetic field.
int GoodCluster(const BaseParticlePropagator &bpg, const reco::PFClusterCollection &pfc, float minep, bool sec=false)
const TrackerGeometry * tracker_
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
tmp
align.sh
Definition: createJobs.py:716
const edm::ESGetToken< GeometricSearchTracker, TrackerRecoGeometryRecord > geomSearchTrackerToken_
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
def move(src, dest)
Definition: eostools.py:511
bool isTrackerPixel(GeomDetEnumerators::SubDetector m)
Definition: Run.h:45
math::XYZTLorentzVector XYZTLorentzVector
Definition: RawParticle.h:25
ib
Definition: cuy.py:661
const XYZTLorentzVector & vertex() const
the vertex fourvector
Definition: RawParticle.h:320
#define LogDebug(id)