CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
MuonSimHitProducer Class Reference
Inheritance diagram for MuonSimHitProducer:
edm::stream::EDProducer<>

Public Member Functions

 MuonSimHitProducer (const edm::ParameterSet &)
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Private Member Functions

void applyMaterialEffects (TrajectoryStateOnSurface &tsosWithdEdx, TrajectoryStateOnSurface &tsos, double radPath, RandomEngineAndDistribution const *, HepPDT::ParticleDataTable const &)
 Simulate material effects in iron (dE/dx, multiple scattering) More...
 
void beginRun (edm::Run const &run, const edm::EventSetup &es) override
 
void produce (edm::Event &, const edm::EventSetup &) override
 
void readParameters (const edm::ParameterSet &, const edm::ParameterSet &, const edm::ParameterSet &)
 

Private Attributes

const CSCGeometrycscGeom
 
const edm::ESGetToken< CSCGeometry, MuonGeometryRecordcscGeometryESToken_
 
bool doGL_
 
bool doL1_
 
bool doL3_
 
const DTGeometrydtGeom
 
const edm::ESGetToken< DTGeometry, MuonGeometryRecorddtGeometryESToken_
 
double fCSC
 
double fDT
 
bool fullPattern_
 
double kCSC
 
double kDT
 
const MagneticFieldmagfield
 
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecordmagneticFieldESToken_
 
const edm::ESGetToken< HepPDT::ParticleDataTable, edm::DefaultRecordparticleDataTableESToken_
 
const PropagatorpropagatorWithMaterial
 
std::unique_ptr< PropagatorpropagatorWithoutMaterial
 
const RPCGeometryrpcGeom
 
const edm::ESGetToken< RPCGeometry, MuonGeometryRecordrpcGeometryESToken_
 
edm::InputTag simMuonLabel
 
edm::EDGetTokenT< std::vector< SimTrack > > simMuonToken
 
edm::InputTag simVertexLabel
 
edm::EDGetTokenT< std::vector< SimVertex > > simVertexToken
 
Chi2MeasurementEstimator theEstimator
 
std::unique_ptr< MaterialEffectstheMaterialEffects
 
MuonServiceProxytheService
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 

Detailed Description

Description: Fast simulation producer of Muon Sim Hits (to be used for realistic Muon reconstruction)

Implementation: <Notes on="" implementation>="">

Definition at line 51 of file MuonSimHitProducer.cc.

Constructor & Destructor Documentation

◆ MuonSimHitProducer()

MuonSimHitProducer::MuonSimHitProducer ( const edm::ParameterSet iConfig)
explicit

Definition at line 110 of file MuonSimHitProducer.cc.

References edm::ParameterSet::getParameter(), MuonServiceProxy_cff::MuonServiceProxy, readParameters(), MuonServiceProxy::Run, simMuonLabel, simMuonToken, simVertexLabel, simVertexToken, and theService.

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>()),
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 }
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_
std::unique_ptr< Propagator > propagatorWithoutMaterial
edm::EDGetTokenT< std::vector< SimVertex > > simVertexToken
const edm::ESGetToken< DTGeometry, MuonGeometryRecord > dtGeometryESToken_
const edm::ESGetToken< RPCGeometry, MuonGeometryRecord > rpcGeometryESToken_
edm::EDGetTokenT< std::vector< SimTrack > > simMuonToken
MuonServiceProxy * theService
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magneticFieldESToken_
void readParameters(const edm::ParameterSet &, const edm::ParameterSet &, const edm::ParameterSet &)
edm::InputTag simMuonLabel
const edm::ESGetToken< CSCGeometry, MuonGeometryRecord > cscGeometryESToken_
edm::InputTag simVertexLabel
Chi2MeasurementEstimator theEstimator

Member Function Documentation

◆ applyMaterialEffects()

void MuonSimHitProducer::applyMaterialEffects ( TrajectoryStateOnSurface tsosWithdEdx,
TrajectoryStateOnSurface tsos,
double  radPath,
RandomEngineAndDistribution const *  random,
HepPDT::ParticleDataTable const &  table 
)
private

Simulate material effects in iron (dE/dx, multiple scattering)

Definition at line 553 of file MuonSimHitProducer.cc.

References fastSimProducer_cff::bremsstrahlung, ALCARECOTkAlJpsiMuMu_cff::charge, TrajectoryStateOnSurface::charge(), fastSimProducer_cff::energyLoss, dqmMemoryStats::float, TrajectoryStateOnSurface::globalMomentum(), TrajectoryStateOnSurface::globalPosition(), PV3DBase< T, PVType, FrameType >::mag2(), magfield, rawparticle::makeMuon(), RawParticle::momentum(), amptDefaultParameters_cff::mu, fastSimProducer_cff::multipleScattering, BaseParticlePropagator::particle(), position, RawParticle::Px(), RawParticle::Py(), RawParticle::Pz(), RawParticle::setMomentum(), RawParticle::setVertex(), mathSSE::sqrt(), TrajectoryStateOnSurface::surface(), TableParser::table, Surface::tangentPlane(), theMaterialEffects, RawParticle::vertex(), PV3DBase< T, PVType, FrameType >::x(), RawParticle::X(), PV3DBase< T, PVType, FrameType >::y(), RawParticle::Y(), PV3DBase< T, PVType, FrameType >::z(), and RawParticle::Z().

Referenced by produce().

557  {
558  // The energy loss simulator
559  EnergyLossSimulator* energyLoss = theMaterialEffects->energyLossSimulator();
560 
561  // The multiple scattering simulator
562  MultipleScatteringSimulator* multipleScattering = theMaterialEffects->multipleScatteringSimulator();
563 
564  // The Muon Bremsstrahlung simulator
565  MuonBremsstrahlungSimulator* bremsstrahlung = theMaterialEffects->muonBremsstrahlungSimulator();
566 
567  // Initialize the Particle position, momentum and energy
568  const Surface& nextSurface = tsos.surface();
569  GlobalPoint gPos = energyLoss ? tsos.globalPosition() : tsosWithdEdx.globalPosition();
570  GlobalVector gMom = energyLoss ? tsos.globalMomentum() : tsosWithdEdx.globalMomentum();
571  double mu = 0.1056583692;
572  double en = std::sqrt(gMom.mag2() + mu * mu);
573 
574  // And now create the Particle
575  XYZTLorentzVector position(gPos.x(), gPos.y(), gPos.z(), 0.);
576  XYZTLorentzVector momentum(gMom.x(), gMom.y(), gMom.z(), en);
577  float charge = (float)(tsos.charge());
578  ParticlePropagator theMuon(rawparticle::makeMuon(charge < 1., momentum, position), nullptr, nullptr, &table);
579 
580  // Recompute the energy loss to get the fluctuations
581  if (energyLoss) {
582  // Difference between with and without dE/dx (average only)
583  // (for corrections once fluctuations are applied)
584  GlobalPoint gPosWithdEdx = tsosWithdEdx.globalPosition();
585  GlobalVector gMomWithdEdx = tsosWithdEdx.globalMomentum();
586  double enWithdEdx = std::sqrt(gMomWithdEdx.mag2() + mu * mu);
587  XYZTLorentzVector deltaPos(
588  gPosWithdEdx.x() - gPos.x(), gPosWithdEdx.y() - gPos.y(), gPosWithdEdx.z() - gPos.z(), 0.);
589  XYZTLorentzVector deltaMom(
590  gMomWithdEdx.x() - gMom.x(), gMomWithdEdx.y() - gMom.y(), gMomWithdEdx.z() - gMom.z(), enWithdEdx - en);
591 
592  // Energy loss in iron (+ fluctuations)
593  energyLoss->updateState(theMuon, radPath, random);
594 
595  // Correcting factors to account for slight differences in material descriptions
596  // (Material description is more accurate in the stepping helix propagator)
597  radPath *= -deltaMom.E() / energyLoss->mostLikelyLoss();
598  double fac = energyLoss->deltaMom().E() / energyLoss->mostLikelyLoss();
599 
600  // Particle momentum & position after energy loss + fluctuation
601  XYZTLorentzVector theNewMomentum = theMuon.particle().momentum() + energyLoss->deltaMom() + fac * deltaMom;
602  XYZTLorentzVector theNewPosition = theMuon.particle().vertex() + fac * deltaPos;
603  fac = (theNewMomentum.E() * theNewMomentum.E() - mu * mu) / theNewMomentum.Vect().Mag2();
604  fac = fac > 0. ? std::sqrt(fac) : 1E-9;
605  theMuon.particle().setMomentum(
606  theNewMomentum.Px() * fac, theNewMomentum.Py() * fac, theNewMomentum.Pz() * fac, theNewMomentum.E());
607  theMuon.particle().setVertex(theNewPosition);
608  }
609 
610  // Does the actual mutliple scattering
611  if (multipleScattering) {
612  // Pass the vector normal to the "next" surface
613  GlobalVector normal = nextSurface.tangentPlane(tsos.globalPosition())->normalVector();
614  multipleScattering->setNormalVector(normal);
615  // Compute the amount of multiple scattering after a given path length
616  multipleScattering->updateState(theMuon, radPath, random);
617  }
618 
619  // Muon Bremsstrahlung
620  if (bremsstrahlung) {
621  // Compute the amount of Muon Bremsstrahlung after given path length
622  bremsstrahlung->updateState(theMuon, radPath, random);
623  }
624 
625  // Fill the propagated state
626  GlobalPoint propagatedPosition(theMuon.particle().X(), theMuon.particle().Y(), theMuon.particle().Z());
627  GlobalVector propagatedMomentum(theMuon.particle().Px(), theMuon.particle().Py(), theMuon.particle().Pz());
628  GlobalTrajectoryParameters propagatedGtp(propagatedPosition, propagatedMomentum, (int)charge, magfield);
629  tsosWithdEdx = TrajectoryStateOnSurface(propagatedGtp, nextSurface);
630 }
const MagneticField * magfield
virtual ConstReferenceCountingPointer< TangentPlane > tangentPlane(const GlobalPoint &) const =0
T z() const
Definition: PV3DBase.h:61
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
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
GlobalPoint globalPosition() const
T sqrt(T t)
Definition: SSEVec.h:19
GlobalVector globalMomentum() const
static int position[264][3]
Definition: ReadPGInfo.cc:289
std::unique_ptr< MaterialEffects > theMaterialEffects
math::XYZTLorentzVector XYZTLorentzVector
Definition: RawParticle.h:25

◆ beginRun()

void MuonSimHitProducer::beginRun ( edm::Run const &  run,
const edm::EventSetup es 
)
overrideprivate

Definition at line 139 of file MuonSimHitProducer.cc.

References Propagator::clone(), cscGeom, cscGeometryESToken_, dtGeom, dtGeometryESToken_, edm::EventSetup::getData(), magfield, magneticFieldESToken_, MuonServiceProxy::propagator(), propagatorWithMaterial, propagatorWithoutMaterial, rpcGeom, rpcGeometryESToken_, SteppingHelixPropagator::setMaterialMode(), theService, and MuonServiceProxy::update().

139  {
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 }
const MagneticField * magfield
const DTGeometry * dtGeom
std::unique_ptr< Propagator > propagatorWithoutMaterial
virtual Propagator * clone() const =0
const edm::ESGetToken< DTGeometry, MuonGeometryRecord > dtGeometryESToken_
const RPCGeometry * rpcGeom
const Propagator * propagatorWithMaterial
const CSCGeometry * cscGeom
bool getData(T &iHolder) const
Definition: EventSetup.h:122
const edm::ESGetToken< RPCGeometry, MuonGeometryRecord > rpcGeometryESToken_
MuonServiceProxy * theService
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > magneticFieldESToken_
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 edm::ESGetToken< CSCGeometry, MuonGeometryRecord > cscGeometryESToken_
void setMaterialMode(bool noMaterial)
Switch for material effects mode: no material effects if true.

◆ produce()

void MuonSimHitProducer::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprivate

Definition at line 163 of file MuonSimHitProducer.cc.

References funct::abs(), alongMomentum, anyDirection, applyMaterialEffects(), PV3DBase< T, PVType, FrameType >::basicVector(), Surface::bounds(), relativeConstraints::chamber, RPCGeometry::chamber(), CSCGeometry::chamber(), DTGeometry::chamber(), DTTopology::channel(), CoreSimTrack::charge(), gather_cfg::cout, GeomDetEnumerators::CSC, cscGeom, MuonServiceProxy::detLayerGeometry(), GeomDetEnumerators::DT, dtGeom, CSCDetId::endcap(), mps_splice::entry, beamvalidation::exit(), JetChargeProducer_cfi::exp, fCSC, fDT, DTTopology::firstChannel(), TrajectoryStateOnSurface::freeTrajectoryState(), GeomDet::geographicalId(), edm::EventSetup::getHandle(), SteppingHelixStateInfo::getStateOnSurface(), TrajectoryStateOnSurface::globalMomentum(), TrajectoryStateOnSurface::globalPosition(), mps_fire::i, triggerObjects_cff::id, RPCGeometry::idToDetUnit(), DTGeometry::idToDetUnit(), CSCGeometry::idToDetUnit(), iEvent, Bounds::inside(), CSCLayerGeometry::inside(), TrajectoryStateOnSurface::isValid(), kCSC, kDT, DTTopology::lastChannel(), DTLayerId::layer(), CSCDetId::layer(), phase1PixelTopology::layer, hgcalTopologyTester_cfi::layers, TrajectoryStateOnSurface::localMomentum(), dqm-mbProfile::log, visualization-live-secondInstance_cfg::m, PV3DBase< T, PVType, FrameType >::mag(), magfield, CoreSimTrack::momentum(), eostools::move(), dqmiodumpmetadata::n, GetRecoTauVFromDQM_MC_cff::next, particleDataTableESToken_, castor_dqm_sourceclient_file_cfg::path, SteppingHelixStateInfo::path(), packedPFCandidateRefMixer_cfi::pf, PV3DBase< T, PVType, FrameType >::phi(), pi, PlaneBuilder::plane(), position, propagatorWithMaterial, propagatorWithoutMaterial, SteppingHelixStateInfo::radPath(), DetId::rawId(), bril_dqm_clientPB-live_cfg::rid, CSCDetId::ring(), relativeConstraints::ring, makeMuonMisalignmentScenario::rot, GeomDetEnumerators::RPCBarrel, GeomDetEnumerators::RPCEndcap, rpcGeom, DTChamberId::sector(), simMuonToken, simVertexToken, HGCalValidator_cfi::simVertices, mathSSE::sqrt(), DTChamberId::station(), relativeConstraints::station, CSCDetId::station(), GeomDet::subDetector(), DTSuperLayerId::superlayer(), GeomDet::surface(), FrontierCondition_GT_autoExpress_cfi::t0, theEstimator, theMaterialEffects, theService, PV3DBase< T, PVType, FrameType >::theta(), Bounds::thickness(), Calorimetry_cff::thickness, GeomDet::toLocal(), SimTrack::trackerSurfaceMomentum(), SimTrack::trackerSurfacePosition(), CoreSimTrack::trackId(), TrajectoryStateOnSurface::transverseCurvature(), CoreSimTrack::type(), Vector3DBase< T, FrameTag >::unit(), SimTrack::vertIndex(), DTChamberId::wheel(), x, HLT_2022v11_cff::xAxis, y, HLT_2022v11_cff::yAxis, z, PV3DBase< T, PVType, FrameType >::z(), and HLT_2022v11_cff::zAxis.

163  {
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  int n = 0;
495  for (std::vector<PSimHit>::const_iterator i = theCSCHits.begin(); i != theCSCHits.end(); i++) {
496  pcsc->push_back(*i);
497  n += 1;
498  }
499  iEvent.put(std::move(pcsc), "MuonCSCHits");
500 
501  std::unique_ptr<edm::PSimHitContainer> pdt(new edm::PSimHitContainer);
502  n = 0;
503  for (std::vector<PSimHit>::const_iterator i = theDTHits.begin(); i != theDTHits.end(); i++) {
504  pdt->push_back(*i);
505  n += 1;
506  }
507  iEvent.put(std::move(pdt), "MuonDTHits");
508 
509  std::unique_ptr<edm::PSimHitContainer> prpc(new edm::PSimHitContainer);
510  n = 0;
511  for (std::vector<PSimHit>::const_iterator i = theRPCHits.begin(); i != theRPCHits.end(); i++) {
512  prpc->push_back(*i);
513  n += 1;
514  }
515  iEvent.put(std::move(prpc), "MuonRPCHits");
516 }
int station() const
Return the station number.
Definition: DTChamberId.h:42
const MagneticField * magfield
const edm::ESGetToken< HepPDT::ParticleDataTable, edm::DefaultRecord > particleDataTableESToken_
TrajectoryStateOnSurface getStateOnSurface(const Surface &surf, bool returnTangentPlane=false) const
const DTGeometry * dtGeom
const math::XYZVectorD & trackerSurfacePosition() const
Definition: SimTrack.h:40
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
T z() const
Definition: PV3DBase.h:61
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
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
int layer() const
Definition: CSCDetId.h:56
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
constexpr std::array< uint8_t, layerIndexSize > layer
const Double_t pi
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
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
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
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:151
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
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:39
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:54
unsigned int trackId() const
Definition: CoreSimTrack.h:31
TkRotation< float > RotationType
int lastChannel() const
Returns the wire number of the last wire.
Definition: DTTopology.h:81
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
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
def move(src, dest)
Definition: eostools.py:511
Global3DVector GlobalVector
Definition: GlobalVector.h:10
Geom::Theta< T > theta() const
Definition: PV3DBase.h:72
const GeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.
Definition: CSCGeometry.cc:89
Chi2MeasurementEstimator theEstimator
def exit(msg="")
const Bounds & bounds() const
Definition: Surface.h:87

◆ readParameters()

void MuonSimHitProducer::readParameters ( const edm::ParameterSet fastMuons,
const edm::ParameterSet fastTracks,
const edm::ParameterSet matEff 
)
private

Definition at line 518 of file MuonSimHitProducer.cc.

References fCSC, fDT, fullPattern_, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), HLT_2022v11_cff::InputTag, kCSC, kDT, simMuonLabel, simVertexLabel, AlCaHLTBitMon_QueryRunRegistry::string, and theMaterialEffects.

Referenced by MuonSimHitProducer().

520  {
521  // Muons
522  std::string _simModuleLabel = fastMuons.getParameter<std::string>("simModuleLabel");
523  std::string _simModuleProcess = fastMuons.getParameter<std::string>("simModuleProcess");
524  simMuonLabel = edm::InputTag(_simModuleLabel, _simModuleProcess);
525  simVertexLabel = edm::InputTag(_simModuleLabel);
526 
527  std::vector<double> simHitIneffDT = fastMuons.getParameter<std::vector<double> >("simHitDTIneffParameters");
528  std::vector<double> simHitIneffCSC = fastMuons.getParameter<std::vector<double> >("simHitCSCIneffParameters");
529  kDT = simHitIneffDT[0];
530  fDT = simHitIneffDT[1];
531  kCSC = simHitIneffCSC[0];
532  fCSC = simHitIneffCSC[1];
533 
534  // Tracks
535  fullPattern_ = fastTracks.getUntrackedParameter<bool>("FullPatternRecognition");
536 
537  // The following should be on LogInfo
538  // std::cout << " MUON SIM HITS: FastSimulation parameters " << std::endl;
539  // std::cout << " ============================================== " << std::endl;
540  // if ( fullPattern_ )
541  // std::cout << " The FULL pattern recognition option is turned ON" << std::endl;
542  // else
543  // std::cout << " The FAST tracking option is turned ON" << std::endl;
544 
545  // Material Effects
546  theMaterialEffects = nullptr;
547  if (matEff.getParameter<bool>("PairProduction") || matEff.getParameter<bool>("Bremsstrahlung") ||
548  matEff.getParameter<bool>("MuonBremsstrahlung") || matEff.getParameter<bool>("EnergyLoss") ||
549  matEff.getParameter<bool>("MultipleScattering"))
550  theMaterialEffects = std::make_unique<MaterialEffects>(matEff);
551 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
T getUntrackedParameter(std::string const &, T const &) const
edm::InputTag simMuonLabel
std::unique_ptr< MaterialEffects > theMaterialEffects
edm::InputTag simVertexLabel

Member Data Documentation

◆ cscGeom

const CSCGeometry* MuonSimHitProducer::cscGeom
private

Definition at line 61 of file MuonSimHitProducer.cc.

Referenced by beginRun(), and produce().

◆ cscGeometryESToken_

const edm::ESGetToken<CSCGeometry, MuonGeometryRecord> MuonSimHitProducer::cscGeometryESToken_
private

Definition at line 99 of file MuonSimHitProducer.cc.

Referenced by beginRun().

◆ doGL_

bool MuonSimHitProducer::doGL_
private

Definition at line 87 of file MuonSimHitProducer.cc.

◆ doL1_

bool MuonSimHitProducer::doL1_
private

Definition at line 87 of file MuonSimHitProducer.cc.

◆ doL3_

bool MuonSimHitProducer::doL3_
private

Definition at line 87 of file MuonSimHitProducer.cc.

◆ dtGeom

const DTGeometry* MuonSimHitProducer::dtGeom
private

Definition at line 60 of file MuonSimHitProducer.cc.

Referenced by beginRun(), and produce().

◆ dtGeometryESToken_

const edm::ESGetToken<DTGeometry, MuonGeometryRecord> MuonSimHitProducer::dtGeometryESToken_
private

Definition at line 98 of file MuonSimHitProducer.cc.

Referenced by beginRun().

◆ fCSC

double MuonSimHitProducer::fCSC
private

Definition at line 76 of file MuonSimHitProducer.cc.

Referenced by produce(), and readParameters().

◆ fDT

double MuonSimHitProducer::fDT
private

Definition at line 74 of file MuonSimHitProducer.cc.

Referenced by produce(), and readParameters().

◆ fullPattern_

bool MuonSimHitProducer::fullPattern_
private

Definition at line 86 of file MuonSimHitProducer.cc.

Referenced by readParameters().

◆ kCSC

double MuonSimHitProducer::kCSC
private

Definition at line 75 of file MuonSimHitProducer.cc.

Referenced by produce(), and readParameters().

◆ kDT

double MuonSimHitProducer::kDT
private

Definition at line 73 of file MuonSimHitProducer.cc.

Referenced by produce(), and readParameters().

◆ magfield

const MagneticField* MuonSimHitProducer::magfield
private

Definition at line 59 of file MuonSimHitProducer.cc.

Referenced by applyMaterialEffects(), beginRun(), and produce().

◆ magneticFieldESToken_

const edm::ESGetToken<MagneticField, IdealMagneticFieldRecord> MuonSimHitProducer::magneticFieldESToken_
private

Definition at line 97 of file MuonSimHitProducer.cc.

Referenced by beginRun().

◆ particleDataTableESToken_

const edm::ESGetToken<HepPDT::ParticleDataTable, edm::DefaultRecord> MuonSimHitProducer::particleDataTableESToken_
private

Definition at line 101 of file MuonSimHitProducer.cc.

Referenced by produce().

◆ propagatorWithMaterial

const Propagator* MuonSimHitProducer::propagatorWithMaterial
private

Definition at line 63 of file MuonSimHitProducer.cc.

Referenced by beginRun(), and produce().

◆ propagatorWithoutMaterial

std::unique_ptr<Propagator> MuonSimHitProducer::propagatorWithoutMaterial
private

Definition at line 64 of file MuonSimHitProducer.cc.

Referenced by beginRun(), and produce().

◆ rpcGeom

const RPCGeometry* MuonSimHitProducer::rpcGeom
private

Definition at line 62 of file MuonSimHitProducer.cc.

Referenced by beginRun(), and produce().

◆ rpcGeometryESToken_

const edm::ESGetToken<RPCGeometry, MuonGeometryRecord> MuonSimHitProducer::rpcGeometryESToken_
private

Definition at line 100 of file MuonSimHitProducer.cc.

Referenced by beginRun().

◆ simMuonLabel

edm::InputTag MuonSimHitProducer::simMuonLabel
private

Definition at line 90 of file MuonSimHitProducer.cc.

Referenced by MuonSimHitProducer(), and readParameters().

◆ simMuonToken

edm::EDGetTokenT<std::vector<SimTrack> > MuonSimHitProducer::simMuonToken
private

Definition at line 94 of file MuonSimHitProducer.cc.

Referenced by MuonSimHitProducer(), and produce().

◆ simVertexLabel

edm::InputTag MuonSimHitProducer::simVertexLabel
private

Definition at line 91 of file MuonSimHitProducer.cc.

Referenced by MuonSimHitProducer(), and readParameters().

◆ simVertexToken

edm::EDGetTokenT<std::vector<SimVertex> > MuonSimHitProducer::simVertexToken
private

Definition at line 95 of file MuonSimHitProducer.cc.

Referenced by MuonSimHitProducer(), and produce().

◆ theEstimator

Chi2MeasurementEstimator MuonSimHitProducer::theEstimator
private

Definition at line 57 of file MuonSimHitProducer.cc.

Referenced by produce().

◆ theMaterialEffects

std::unique_ptr<MaterialEffects> MuonSimHitProducer::theMaterialEffects
private

Definition at line 66 of file MuonSimHitProducer.cc.

Referenced by applyMaterialEffects(), produce(), and readParameters().

◆ theService

MuonServiceProxy* MuonSimHitProducer::theService
private

Definition at line 56 of file MuonSimHitProducer.cc.

Referenced by beginRun(), MuonSimHitProducer(), and produce().