CMS 3D CMS Logo

MuonSimHitProducer.cc
Go to the documentation of this file.
1 //
2 // Package: MuonSimHitProducer
3 // Class: MuonSimHitProducer
4 //
14 //
15 // Original Author: Martijn Mulders/Matthew Jones
16 // Created: Wed Jul 30 11:37:24 CET 2007
17 // Working: Fri Nov 9 09:39:33 CST 2007
18 //
19 // $Id: MuonSimHitProducer.cc,v 1.36 2011/10/07 08:25:42 aperrott Exp $
20 //
21 //
22 
50 
52 public:
53  explicit MuonSimHitProducer(const edm::ParameterSet&);
54 
55 private:
58 
64  std::unique_ptr<Propagator> propagatorWithoutMaterial;
65 
66  std::unique_ptr<MaterialEffects> theMaterialEffects;
67 
68  void beginRun(edm::Run const& run, const edm::EventSetup& es) override;
69  void produce(edm::Event&, const edm::EventSetup&) override;
71 
72  // Parameters to emulate the muonSimHit association inefficiency due to delta's
73  double kDT;
74  double fDT;
75  double kCSC;
76  double fCSC;
77 
81  double radPath,
84 
85  // ----------- parameters ----------------------------
87  bool doL1_, doL3_, doGL_;
88 
89  // tags
92 
93  // tokens
96 
102 };
103 
104 //for debug only
105 //#define FAMOS_DEBUG
106 
107 //
108 // constructors and destructor
109 //
111  : theEstimator(iConfig.getParameter<double>("Chi2EstimatorCut")),
112  propagatorWithoutMaterial(nullptr),
113  magneticFieldESToken_(esConsumes<edm::Transition::BeginRun>()),
114  dtGeometryESToken_(esConsumes<edm::Transition::BeginRun>(edm::ESInputTag("", "MisAligned"))),
115  cscGeometryESToken_(esConsumes<edm::Transition::BeginRun>(edm::ESInputTag("", "MisAligned"))),
116  rpcGeometryESToken_(esConsumes<edm::Transition::BeginRun>()),
117  particleDataTableESToken_(esConsumes()) {
118  // Read relevant parameters
120  iConfig.getParameter<edm::ParameterSet>("TRACKS"),
121  iConfig.getParameter<edm::ParameterSet>("MaterialEffectsForMuons"));
122 
123  //
124  // register your products ... need to declare at least one possible product...
125  //
126  produces<edm::PSimHitContainer>("MuonCSCHits");
127  produces<edm::PSimHitContainer>("MuonDTHits");
128  produces<edm::PSimHitContainer>("MuonRPCHits");
129 
130  edm::ParameterSet serviceParameters = iConfig.getParameter<edm::ParameterSet>("ServiceParameters");
131  theService = new MuonServiceProxy(serviceParameters, consumesCollector(), MuonServiceProxy::UseEventSetupIn::Run);
132 
133  // consumes
134  simMuonToken = consumes<std::vector<SimTrack> >(simMuonLabel);
135  simVertexToken = consumes<std::vector<SimVertex> >(simVertexLabel);
136 }
137 
138 // ---- method called once each job just before starting event loop ----
140  //services
145 
146  bool duringEvent = false;
147  theService->update(es, duringEvent);
148 
149  // A few propagators
150  propagatorWithMaterial = &(*(theService->propagator("SteppingHelixPropagatorAny")));
152  SteppingHelixPropagator* SHpropagator =
153  dynamic_cast<SteppingHelixPropagator*>(propagatorWithoutMaterial.get()); // Beuark!
154  SHpropagator->setMaterialMode(true); // switches OFF material effects;
155 }
156 
157 //
158 // member functions
159 //
160 
161 // ------------ method called to produce the data ------------
162 
164  auto pdg = iSetup.getHandle(particleDataTableESToken_);
165 
166  RandomEngineAndDistribution random(iEvent.streamID());
167 
168  MuonPatternRecoDumper dumper;
169 
172  std::vector<PSimHit> theCSCHits;
173  std::vector<PSimHit> theDTHits;
174  std::vector<PSimHit> theRPCHits;
175 
177  iEvent.getByToken(simMuonToken, simMuons);
178  iEvent.getByToken(simVertexToken, simVertices);
179 
180  for (unsigned int itrk = 0; itrk < simMuons->size(); itrk++) {
181  const SimTrack& mySimTrack = (*simMuons)[itrk];
182  math::XYZTLorentzVector mySimP4(
183  mySimTrack.momentum().x(), mySimTrack.momentum().y(), mySimTrack.momentum().z(), mySimTrack.momentum().t());
184 
185  // Decaying hadrons are now in the list, and so are their muon daughter
186  // Ignore the hadrons here.
187  int pid = mySimTrack.type();
188  if (abs(pid) != 13 && abs(pid) != 1000024)
189  continue;
190  double t0 = 0;
191  GlobalPoint initialPosition;
192  int ivert = mySimTrack.vertIndex();
193  if (ivert >= 0) {
194  t0 = (*simVertices)[ivert].position().t();
195  GlobalPoint xyzzy((*simVertices)[ivert].position().x(),
196  (*simVertices)[ivert].position().y(),
197  (*simVertices)[ivert].position().z());
198  initialPosition = xyzzy;
199  }
200  //
201  // Presumably t0 has dimensions of cm if not mm?
202  // Convert to ns for internal calculations.
203  // I wonder where we should get c from?
204  //
205  double tof = t0 / 29.98;
206 
207 #ifdef FAMOS_DEBUG
208  std::cout << " ===> MuonSimHitProducer::reconstruct() found SIMTRACK - pid = " << pid;
209  std::cout << " : pT = " << mySimP4.Pt() << ", eta = " << mySimP4.Eta() << ", phi = " << mySimP4.Phi() << std::endl;
210 #endif
211 
212  //
213  // Produce muons sim hits starting from undecayed simulated muons
214  //
215 
216  GlobalPoint startingPosition(mySimTrack.trackerSurfacePosition().x(),
217  mySimTrack.trackerSurfacePosition().y(),
218  mySimTrack.trackerSurfacePosition().z());
219  GlobalVector startingMomentum(mySimTrack.trackerSurfaceMomentum().x(),
220  mySimTrack.trackerSurfaceMomentum().y(),
221  mySimTrack.trackerSurfaceMomentum().z());
222  //
223  // Crap... there's no time-of-flight to the trackerSurfacePosition()...
224  // So, this will be wrong when the curvature can't be neglected, but that
225  // will be rather seldom... May as well ignore the mass too.
226  //
227  GlobalVector dtracker = startingPosition - initialPosition;
228  tof += dtracker.mag() / 29.98;
229 
230 #ifdef FAMOS_DEBUG
231  std::cout << " the Muon START position " << startingPosition << std::endl;
232  std::cout << " the Muon START momentum " << startingMomentum << std::endl;
233 #endif
234 
235  //
236  // Some magic to define a TrajectoryStateOnSurface
237  //
238  PlaneBuilder pb;
239  GlobalVector zAxis = startingMomentum.unit();
240  GlobalVector yAxis(zAxis.y(), -zAxis.x(), 0);
241  GlobalVector xAxis = yAxis.cross(zAxis);
243  PlaneBuilder::ReturnType startingPlane = pb.plane(startingPosition, rot);
244  GlobalTrajectoryParameters gtp(startingPosition, startingMomentum, (int)mySimTrack.charge(), magfield);
245  TrajectoryStateOnSurface startingState(gtp, *startingPlane);
246 
247  std::vector<const DetLayer*> navLayers;
248  if (fabs(startingState.globalMomentum().eta()) > 4.5) {
249  navLayers = navigation.compatibleEndcapLayers(*(startingState.freeState()), alongMomentum);
250  } else {
251  navLayers = navigation.compatibleLayers(*(startingState.freeState()), alongMomentum);
252  }
253  /*
254  edm::ESHandle<Propagator> propagator =
255  theService->propagator("SteppingHelixPropagatorAny");
256  */
257 
258  if (navLayers.empty())
259  continue;
260 
261 #ifdef FAMOS_DEBUG
262  std::cout << "Found " << navLayers.size() << " compatible DetLayers..." << std::endl;
263 #endif
264 
265  TrajectoryStateOnSurface propagatedState = startingState;
266  for (unsigned int ilayer = 0; ilayer < navLayers.size(); ilayer++) {
267 #ifdef FAMOS_DEBUG
268  std::cout << "Propagating to layer " << ilayer << " " << dumper.dumpLayer(navLayers[ilayer]) << std::endl;
269 #endif
270 
271  std::vector<DetWithState> comps =
272  navLayers[ilayer]->compatibleDets(propagatedState, *propagatorWithMaterial, theEstimator);
273  if (comps.empty())
274  continue;
275 
276 #ifdef FAMOS_DEBUG
277  std::cout << "Propagating " << propagatedState << std::endl;
278 #endif
279 
280  // Starting momentum
281  double pi = propagatedState.globalMomentum().mag();
282 
283  // Propagate with material effects (dE/dx average only)
284  SteppingHelixStateInfo shsStart(*(propagatedState.freeTrajectoryState()));
285  SteppingHelixStateInfo shsDest;
287  ->propagate(shsStart, navLayers[ilayer]->surface(), shsDest);
288  std::pair<TrajectoryStateOnSurface, double> next(shsDest.getStateOnSurface(navLayers[ilayer]->surface()),
289  shsDest.path());
290  // No need to continue if there is no valid propagation available.
291  // This happens rarely (~0.1% of ttbar events)
292  if (!next.first.isValid())
293  continue;
294  // This is the estimate of the number of radiation lengths traversed,
295  // together with the total path length
296  double radPath = shsDest.radPath();
297  double pathLength = next.second;
298 
299  // Now propagate without dE/dx (average)
300  // [To add the dE/dx fluctuations to the actual dE/dx]
301  std::pair<TrajectoryStateOnSurface, double> nextNoMaterial =
302  propagatorWithoutMaterial->propagateWithPath(propagatedState, navLayers[ilayer]->surface());
303 
304  // Update the propagated state
305  propagatedState = next.first;
306  double pf = propagatedState.globalMomentum().mag();
307 
308  // Insert dE/dx fluctuations and multiple scattering
309  // Skip this step if nextNoMaterial.first is not valid
310  // This happens rarely (~0.02% of ttbar events)
311  if (theMaterialEffects && nextNoMaterial.first.isValid())
312  applyMaterialEffects(propagatedState, nextNoMaterial.first, radPath, &random, *pdg);
313  // Check that the 'shaken' propagatedState is still valid, otherwise continue
314  if (!propagatedState.isValid())
315  continue;
316  // (No evidence that this ever happens)
317  //
318  // Consider this... 1 GeV muon has a velocity that is only 0.5% slower than c...
319  // We probably can safely ignore the mass for anything that makes it out to the
320  // muon chambers.
321  //
322  double pavg = 0.5 * (pi + pf);
323  double m = mySimP4.M();
324  double rbeta = sqrt(1 + m * m / (pavg * pavg)) / 29.98;
325  double dtof = pathLength * rbeta;
326 
327 #ifdef FAMOS_DEBUG
328  std::cout << "Propagated to next surface... path length = " << pathLength << " cm, dTOF = " << dtof << " ns"
329  << std::endl;
330 #endif
331 
332  tof += dtof;
333 
334  for (unsigned int icomp = 0; icomp < comps.size(); icomp++) {
335  const GeomDet* gd = comps[icomp].first;
336  if (gd->subDetector() == GeomDetEnumerators::DT) {
338  const DTChamber* chamber = dtGeom->chamber(id);
339  std::vector<const DTSuperLayer*> superlayer = chamber->superLayers();
340  for (unsigned int isl = 0; isl < superlayer.size(); isl++) {
341  std::vector<const DTLayer*> layer = superlayer[isl]->layers();
342  for (unsigned int ilayer = 0; ilayer < layer.size(); ilayer++) {
343  DTLayerId lid = layer[ilayer]->id();
344 #ifdef FAMOS_DEBUG
345  std::cout << " Extrapolated to DT (" << lid.wheel() << "," << lid.station() << "," << lid.sector()
346  << "," << lid.superlayer() << "," << lid.layer() << ")" << std::endl;
347 #endif
348 
349  const GeomDetUnit* det = dtGeom->idToDetUnit(lid);
350 
351  HelixArbitraryPlaneCrossing crossing(propagatedState.globalPosition().basicVector(),
352  propagatedState.globalMomentum().basicVector(),
353  propagatedState.transverseCurvature(),
354  anyDirection);
355  std::pair<bool, double> path = crossing.pathLength(det->surface());
356  if (!path.first)
357  continue;
358  LocalPoint lpos = det->toLocal(GlobalPoint(crossing.position(path.second)));
359  if (!det->surface().bounds().inside(lpos))
360  continue;
361  const DTTopology& dtTopo = layer[ilayer]->specificTopology();
362  int wire = dtTopo.channel(lpos);
363  if (wire - dtTopo.firstChannel() < 0 || wire - dtTopo.lastChannel() > 0)
364  continue;
365  // no drift cell here (on the chamber edge or just outside)
366  // this hit would otherwise be discarded downstream in the digitizer
367 
368  DTWireId wid(lid, wire);
369  double thickness = det->surface().bounds().thickness();
370  LocalVector lmom = det->toLocal(GlobalVector(crossing.direction(path.second)));
371  lmom = lmom.unit() * propagatedState.localMomentum().mag();
372 
373  // Factor that takes into account the (rec)hits lost because of delta's, etc.:
374  // (Not fully satisfactory patch, but it seems to work...)
375  double pmu = lmom.mag();
376  double theDTHitIneff = pmu > 0 ? exp(kDT * log(pmu) + fDT) : 0.;
377  if (random.flatShoot() < theDTHitIneff)
378  continue;
379 
380  double eloss = 0;
381  double pz = fabs(lmom.z());
382  LocalPoint entry = lpos - 0.5 * thickness * lmom / pz;
383  LocalPoint exit = lpos + 0.5 * thickness * lmom / pz;
384  double dtof = path.second * rbeta;
385  int trkid = mySimTrack.trackId();
386  unsigned int id = wid.rawId();
387  short unsigned int processType = 2;
388  PSimHit hit(
389  entry, exit, lmom.mag(), tof + dtof, eloss, pid, id, trkid, lmom.theta(), lmom.phi(), processType);
390  theDTHits.push_back(hit);
391  }
392  }
393  } else if (gd->subDetector() == GeomDetEnumerators::CSC) {
394  CSCDetId id(gd->geographicalId());
395  const CSCChamber* chamber = cscGeom->chamber(id);
396  std::vector<const CSCLayer*> layers = chamber->layers();
397  for (unsigned int ilayer = 0; ilayer < layers.size(); ilayer++) {
398  CSCDetId lid = layers[ilayer]->id();
399 
400 #ifdef FAMOS_DEBUG
401  std::cout << " Extrapolated to CSC (" << lid.endcap() << "," << lid.ring() << "," << lid.station() << ","
402  << lid.layer() << ")" << std::endl;
403 #endif
404 
405  const GeomDetUnit* det = cscGeom->idToDetUnit(lid);
406  HelixArbitraryPlaneCrossing crossing(propagatedState.globalPosition().basicVector(),
407  propagatedState.globalMomentum().basicVector(),
408  propagatedState.transverseCurvature(),
409  anyDirection);
410  std::pair<bool, double> path = crossing.pathLength(det->surface());
411  if (!path.first)
412  continue;
413  LocalPoint lpos = det->toLocal(GlobalPoint(crossing.position(path.second)));
414  // For CSCs the Bounds are for chamber frames not gas regions
415  // if ( ! det->surface().bounds().inside(lpos) ) continue;
416  // New function knows where the 'active' volume is:
417  const CSCLayerGeometry* laygeom = layers[ilayer]->geometry();
418  if (!laygeom->inside(lpos))
419  continue;
420  //double thickness = laygeom->thickness(); gives number which is about 20 times too big
421  double thickness = det->surface().bounds().thickness(); // this one works much better...
422  LocalVector lmom = det->toLocal(GlobalVector(crossing.direction(path.second)));
423  lmom = lmom.unit() * propagatedState.localMomentum().mag();
424 
425  // Factor that takes into account the (rec)hits lost because of delta's, etc.:
426  // (Not fully satisfactory patch, but it seems to work...)
427  double pmu = lmom.mag();
428  double theCSCHitIneff = pmu > 0 ? exp(kCSC * log(pmu) + fCSC) : 0.;
429  // Take into account the different geometry in ME11:
430  if (id.station() == 1 && id.ring() == 1)
431  theCSCHitIneff = theCSCHitIneff * 0.442;
432  if (random.flatShoot() < theCSCHitIneff)
433  continue;
434 
435  double eloss = 0;
436  double pz = fabs(lmom.z());
437  LocalPoint entry = lpos - 0.5 * thickness * lmom / pz;
438  LocalPoint exit = lpos + 0.5 * thickness * lmom / pz;
439  double dtof = path.second * rbeta;
440  int trkid = mySimTrack.trackId();
441  unsigned int id = lid.rawId();
442  short unsigned int processType = 2;
443  PSimHit hit(
444  entry, exit, lmom.mag(), tof + dtof, eloss, pid, id, trkid, lmom.theta(), lmom.phi(), processType);
445  theCSCHits.push_back(hit);
446  }
447  } else if (gd->subDetector() == GeomDetEnumerators::RPCBarrel ||
449  RPCDetId id(gd->geographicalId());
450  const RPCChamber* chamber = rpcGeom->chamber(id);
451  std::vector<const RPCRoll*> roll = chamber->rolls();
452  for (unsigned int iroll = 0; iroll < roll.size(); iroll++) {
453  RPCDetId rid = roll[iroll]->id();
454 
455 #ifdef FAMOS_DEBUG
456  std::cout << " Extrapolated to RPC (" << rid.ring() << "," << rid.station() << "," << rid.sector() << ","
457  << rid.subsector() << "," << rid.layer() << "," << rid.roll() << ")" << std::endl;
458 #endif
459 
460  const GeomDetUnit* det = rpcGeom->idToDetUnit(rid);
461  HelixArbitraryPlaneCrossing crossing(propagatedState.globalPosition().basicVector(),
462  propagatedState.globalMomentum().basicVector(),
463  propagatedState.transverseCurvature(),
464  anyDirection);
465  std::pair<bool, double> path = crossing.pathLength(det->surface());
466  if (!path.first)
467  continue;
468  LocalPoint lpos = det->toLocal(GlobalPoint(crossing.position(path.second)));
469  if (!det->surface().bounds().inside(lpos))
470  continue;
471  double thickness = det->surface().bounds().thickness();
472  LocalVector lmom = det->toLocal(GlobalVector(crossing.direction(path.second)));
473  lmom = lmom.unit() * propagatedState.localMomentum().mag();
474  double eloss = 0;
475  double pz = fabs(lmom.z());
476  LocalPoint entry = lpos - 0.5 * thickness * lmom / pz;
477  LocalPoint exit = lpos + 0.5 * thickness * lmom / pz;
478  double dtof = path.second * rbeta;
479  int trkid = mySimTrack.trackId();
480  unsigned int id = rid.rawId();
481  short unsigned int processType = 2;
482  PSimHit hit(
483  entry, exit, lmom.mag(), tof + dtof, eloss, pid, id, trkid, lmom.theta(), lmom.phi(), processType);
484  theRPCHits.push_back(hit);
485  }
486  } else {
487  std::cout << "Extrapolated to unknown subdetector '" << gd->subDetector() << "'..." << std::endl;
488  }
489  }
490  }
491  }
492 
493  std::unique_ptr<edm::PSimHitContainer> pcsc(new edm::PSimHitContainer);
494  for (std::vector<PSimHit>::const_iterator i = theCSCHits.begin(); i != theCSCHits.end(); i++) {
495  pcsc->push_back(*i);
496  }
497  iEvent.put(std::move(pcsc), "MuonCSCHits");
498 
499  std::unique_ptr<edm::PSimHitContainer> pdt(new edm::PSimHitContainer);
500  for (std::vector<PSimHit>::const_iterator i = theDTHits.begin(); i != theDTHits.end(); i++) {
501  pdt->push_back(*i);
502  }
503  iEvent.put(std::move(pdt), "MuonDTHits");
504 
505  std::unique_ptr<edm::PSimHitContainer> prpc(new edm::PSimHitContainer);
506  for (std::vector<PSimHit>::const_iterator i = theRPCHits.begin(); i != theRPCHits.end(); i++) {
507  prpc->push_back(*i);
508  }
509  iEvent.put(std::move(prpc), "MuonRPCHits");
510 }
511 
513  const edm::ParameterSet& fastTracks,
514  const edm::ParameterSet& matEff) {
515  // Muons
516  std::string _simModuleLabel = fastMuons.getParameter<std::string>("simModuleLabel");
517  std::string _simModuleProcess = fastMuons.getParameter<std::string>("simModuleProcess");
518  simMuonLabel = edm::InputTag(_simModuleLabel, _simModuleProcess);
519  simVertexLabel = edm::InputTag(_simModuleLabel);
520 
521  std::vector<double> simHitIneffDT = fastMuons.getParameter<std::vector<double> >("simHitDTIneffParameters");
522  std::vector<double> simHitIneffCSC = fastMuons.getParameter<std::vector<double> >("simHitCSCIneffParameters");
523  kDT = simHitIneffDT[0];
524  fDT = simHitIneffDT[1];
525  kCSC = simHitIneffCSC[0];
526  fCSC = simHitIneffCSC[1];
527 
528  // Tracks
529  fullPattern_ = fastTracks.getUntrackedParameter<bool>("FullPatternRecognition");
530 
531  // The following should be on LogInfo
532  // std::cout << " MUON SIM HITS: FastSimulation parameters " << std::endl;
533  // std::cout << " ============================================== " << std::endl;
534  // if ( fullPattern_ )
535  // std::cout << " The FULL pattern recognition option is turned ON" << std::endl;
536  // else
537  // std::cout << " The FAST tracking option is turned ON" << std::endl;
538 
539  // Material Effects
540  theMaterialEffects = nullptr;
541  if (matEff.getParameter<bool>("PairProduction") || matEff.getParameter<bool>("Bremsstrahlung") ||
542  matEff.getParameter<bool>("MuonBremsstrahlung") || matEff.getParameter<bool>("EnergyLoss") ||
543  matEff.getParameter<bool>("MultipleScattering"))
544  theMaterialEffects = std::make_unique<MaterialEffects>(matEff);
545 }
546 
549  double radPath,
550  RandomEngineAndDistribution const* random,
552  // The energy loss simulator
553  EnergyLossSimulator* energyLoss = theMaterialEffects->energyLossSimulator();
554 
555  // The multiple scattering simulator
556  MultipleScatteringSimulator* multipleScattering = theMaterialEffects->multipleScatteringSimulator();
557 
558  // The Muon Bremsstrahlung simulator
559  MuonBremsstrahlungSimulator* bremsstrahlung = theMaterialEffects->muonBremsstrahlungSimulator();
560 
561  // Initialize the Particle position, momentum and energy
562  const Surface& nextSurface = tsos.surface();
563  GlobalPoint gPos = energyLoss ? tsos.globalPosition() : tsosWithdEdx.globalPosition();
564  GlobalVector gMom = energyLoss ? tsos.globalMomentum() : tsosWithdEdx.globalMomentum();
565  double mu = 0.1056583692;
566  double en = std::sqrt(gMom.mag2() + mu * mu);
567 
568  // And now create the Particle
569  XYZTLorentzVector position(gPos.x(), gPos.y(), gPos.z(), 0.);
570  XYZTLorentzVector momentum(gMom.x(), gMom.y(), gMom.z(), en);
571  float charge = (float)(tsos.charge());
572  ParticlePropagator theMuon(rawparticle::makeMuon(charge < 1., momentum, position), nullptr, nullptr, &table);
573 
574  // Recompute the energy loss to get the fluctuations
575  if (energyLoss) {
576  // Difference between with and without dE/dx (average only)
577  // (for corrections once fluctuations are applied)
578  GlobalPoint gPosWithdEdx = tsosWithdEdx.globalPosition();
579  GlobalVector gMomWithdEdx = tsosWithdEdx.globalMomentum();
580  double enWithdEdx = std::sqrt(gMomWithdEdx.mag2() + mu * mu);
581  XYZTLorentzVector deltaPos(
582  gPosWithdEdx.x() - gPos.x(), gPosWithdEdx.y() - gPos.y(), gPosWithdEdx.z() - gPos.z(), 0.);
583  XYZTLorentzVector deltaMom(
584  gMomWithdEdx.x() - gMom.x(), gMomWithdEdx.y() - gMom.y(), gMomWithdEdx.z() - gMom.z(), enWithdEdx - en);
585 
586  // Energy loss in iron (+ fluctuations)
587  energyLoss->updateState(theMuon, radPath, random);
588 
589  // Correcting factors to account for slight differences in material descriptions
590  // (Material description is more accurate in the stepping helix propagator)
591  radPath *= -deltaMom.E() / energyLoss->mostLikelyLoss();
592  double fac = energyLoss->deltaMom().E() / energyLoss->mostLikelyLoss();
593 
594  // Particle momentum & position after energy loss + fluctuation
595  XYZTLorentzVector theNewMomentum = theMuon.particle().momentum() + energyLoss->deltaMom() + fac * deltaMom;
596  XYZTLorentzVector theNewPosition = theMuon.particle().vertex() + fac * deltaPos;
597  fac = (theNewMomentum.E() * theNewMomentum.E() - mu * mu) / theNewMomentum.Vect().Mag2();
598  fac = fac > 0. ? std::sqrt(fac) : 1E-9;
599  theMuon.particle().setMomentum(
600  theNewMomentum.Px() * fac, theNewMomentum.Py() * fac, theNewMomentum.Pz() * fac, theNewMomentum.E());
601  theMuon.particle().setVertex(theNewPosition);
602  }
603 
604  // Does the actual mutliple scattering
605  if (multipleScattering) {
606  // Pass the vector normal to the "next" surface
607  GlobalVector normal = nextSurface.tangentPlane(tsos.globalPosition())->normalVector();
608  multipleScattering->setNormalVector(normal);
609  // Compute the amount of multiple scattering after a given path length
610  multipleScattering->updateState(theMuon, radPath, random);
611  }
612 
613  // Muon Bremsstrahlung
614  if (bremsstrahlung) {
615  // Compute the amount of Muon Bremsstrahlung after given path length
616  bremsstrahlung->updateState(theMuon, radPath, random);
617  }
618 
619  // Fill the propagated state
620  GlobalPoint propagatedPosition(theMuon.particle().X(), theMuon.particle().Y(), theMuon.particle().Z());
621  GlobalVector propagatedMomentum(theMuon.particle().Px(), theMuon.particle().Py(), theMuon.particle().Pz());
622  GlobalTrajectoryParameters propagatedGtp(propagatedPosition, propagatedMomentum, (int)charge, magfield);
623  tsosWithdEdx = TrajectoryStateOnSurface(propagatedGtp, nextSurface);
624 }
625 
626 //define this as a plug-in
MuonSimHitProducer(const edm::ParameterSet &)
void setMomentum(const XYZTLorentzVector &vtx)
set the momentum
Definition: RawParticle.h:328
int station() const
Return the station number.
Definition: DTChamberId.h:42
const MagneticField * magfield
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
const edm::ESGetToken< HepPDT::ParticleDataTable, edm::DefaultRecord > particleDataTableESToken_
TrajectoryStateOnSurface getStateOnSurface(const Surface &surf, bool returnTangentPlane=false) const
const DTGeometry * dtGeom
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
double Pz() const
z of the momentum
Definition: RawParticle.h:303
virtual ConstReferenceCountingPointer< TangentPlane > tangentPlane(const GlobalPoint &) const =0
const math::XYZVectorD & trackerSurfacePosition() const
Definition: SimTrack.h:40
static const char * normal
Definition: DMRtrends.cc:35
std::unique_ptr< Propagator > propagatorWithoutMaterial
const CSCChamber * chamber(CSCDetId id) const
Return the chamber corresponding to given DetId.
Definition: CSCGeometry.cc:100
virtual bool inside(const Local3DPoint &) const =0
Determine if the point is inside the bounds.
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
Definition: GeomDet.h:58
HepPDT::ParticleDataTable ParticleDataTable
T z() const
Definition: PV3DBase.h:61
double Z() const
z of vertex
Definition: RawParticle.h:288
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
edm::ESHandle< MuonDetLayerGeometry > detLayerGeometry() const
get the detLayer geometry
ReturnType plane(const PositionType &pos, const RotationType &rot) const
Definition: PlaneBuilder.h:21
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
virtual Propagator * clone() const =0
const GeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.
Definition: RPCGeometry.cc:30
edm::EDGetTokenT< std::vector< SimVertex > > simVertexToken
const XYZTLorentzVector & momentum() const
the momentum fourvector
Definition: RawParticle.h:321
int layer() const
Definition: CSCDetId.h:56
RawParticle makeMuon(bool isParticle, const math::XYZTLorentzVector &p, const math::XYZTLorentzVector &xStart)
Definition: makeMuon.cc:20
T mag2() const
Definition: PV3DBase.h:63
const SurfaceType & surface() const
void applyMaterialEffects(TrajectoryStateOnSurface &tsosWithdEdx, TrajectoryStateOnSurface &tsos, double radPath, RandomEngineAndDistribution const *, HepPDT::ParticleDataTable const &)
Simulate material effects in iron (dE/dx, multiple scattering)
float charge() const
charge
Definition: CoreSimTrack.cc:17
const math::XYZTLorentzVectorD & momentum() const
Definition: CoreSimTrack.h:19
const edm::ESGetToken< DTGeometry, MuonGeometryRecord > dtGeometryESToken_
const Double_t pi
T getUntrackedParameter(std::string const &, T const &) const
const RPCGeometry * rpcGeom
int firstChannel() const
Returns the wire number of the first wire.
Definition: DTTopology.h:79
virtual float thickness() const =0
int type() const
particle type (HEP PDT convension)
Definition: CoreSimTrack.h:22
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
int vertIndex() const
index of the vertex in the Event container (-1 if no vertex)
Definition: SimTrack.h:33
int iEvent
Definition: GenABIO.cc:224
GlobalPoint globalPosition() const
const Propagator * propagatorWithMaterial
T sqrt(T t)
Definition: SSEVec.h:19
double Py() const
y of the momentum
Definition: RawParticle.h:300
FreeTrajectoryState const * freeTrajectoryState(bool withErrors=true) const
T mag() const
Definition: PV3DBase.h:64
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const CSCGeometry * cscGeom
Transition
Definition: Transition.h:12
RawParticle const & particle() const
The particle being propagated.
void beginRun(edm::Run const &run, const edm::EventSetup &es) override
void produce(edm::Event &, const edm::EventSetup &) override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const GeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.
Definition: DTGeometry.cc:75
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:64
const BasicVectorType & basicVector() const
Definition: PV3DBase.h:53
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
const edm::ESGetToken< RPCGeometry, MuonGeometryRecord > rpcGeometryESToken_
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
int superlayer() const
Return the superlayer number (deprecated method name)
int station() const
Definition: CSCDetId.h:79
const RPCChamber * chamber(RPCDetId id) const
Definition: RPCGeometry.cc:46
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
int endcap() const
Definition: CSCDetId.h:85
virtual SubDetector subDetector() const
Which subdetector.
Definition: GeomDet.cc:38
GlobalVector globalMomentum() const
edm::EDGetTokenT< std::vector< SimTrack > > simMuonToken
MuonServiceProxy * theService
int layer() const
Return the layer number.
Definition: DTLayerId.h:42
double Y() const
y of vertex
Definition: RawParticle.h:287
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magneticFieldESToken_
void readParameters(const edm::ParameterSet &, const edm::ParameterSet &, const edm::ParameterSet &)
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:39
HLT enums.
int sector() const
Definition: DTChamberId.h:49
static int position[264][3]
Definition: ReadPGInfo.cc:289
int channel(const LocalPoint &p) const override
Definition: DTTopology.cc:58
unsigned int trackId() const
Definition: CoreSimTrack.h:31
TkRotation< float > RotationType
edm::InputTag simMuonLabel
double X() const
x of vertex
Definition: RawParticle.h:286
int lastChannel() const
Returns the wire number of the last wire.
Definition: DTTopology.h:81
edm::ESHandle< Propagator > propagator(std::string propagatorName) const
get the propagator
void update(const edm::EventSetup &setup, bool duringEvent=true)
update the services each event
const math::XYZTLorentzVectorD & trackerSurfaceMomentum() const
Definition: SimTrack.h:42
std::vector< PSimHit > PSimHitContainer
Vector3DBase unit() const
Definition: Vector3DBase.h:54
int ring() const
Definition: CSCDetId.h:68
const edm::ESGetToken< CSCGeometry, MuonGeometryRecord > cscGeometryESToken_
bool inside(const Local3DPoint &, const LocalError &, float scale=1.f) const override
std::unique_ptr< MaterialEffects > theMaterialEffects
const DTChamber * chamber(const DTChamberId &id) const
Return a DTChamber given its id.
Definition: DTGeometry.cc:90
void setMaterialMode(bool noMaterial)
Switch for material effects mode: no material effects if true.
void setVertex(const XYZTLorentzVector &vtx)
set the vertex
Definition: RawParticle.h:325
def move(src, dest)
Definition: eostools.py:511
Definition: Run.h:45
Global3DVector GlobalVector
Definition: GlobalVector.h:10
Geom::Theta< T > theta() const
Definition: PV3DBase.h:72
math::XYZTLorentzVector XYZTLorentzVector
Definition: RawParticle.h:25
const XYZTLorentzVector & vertex() const
the vertex fourvector
Definition: RawParticle.h:320
double Px() const
x of the momentum
Definition: RawParticle.h:297
const GeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.
Definition: CSCGeometry.cc:89
edm::InputTag simVertexLabel
Chi2MeasurementEstimator theEstimator
def exit(msg="")
const Bounds & bounds() const
Definition: Surface.h:87