CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
MuonSimHitProducer Class Reference

#include <FastSimulation/MuonSimHitProducer/src/MuonSimHitProducer.cc>

Inheritance diagram for MuonSimHitProducer:
edm::EDProducer edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

 MuonSimHitProducer (const edm::ParameterSet &)
 
 ~MuonSimHitProducer ()
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
virtual ~EDProducer ()
 
- Public Member Functions inherited from edm::ProducerBase
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription
const &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 
- Public Member Functions inherited from edm::EDConsumerBase
 EDConsumerBase ()
 
ProductHolderIndex indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndex > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndex > &) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Member Functions

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

Private Attributes

const CSCGeometrycscGeom
 
bool doGL_
 
bool doL1_
 
bool doL3_
 
const DTGeometrydtGeom
 
double fCSC
 
double fDT
 
bool fullPattern_
 
double kCSC
 
double kDT
 
const MagneticFieldmagfield
 
const PropagatorpropagatorWithMaterial
 
PropagatorpropagatorWithoutMaterial
 
const RandomEnginerandom
 
const RPCGeometryrpcGeom
 
Chi2MeasurementEstimator theEstimator
 
MaterialEffectstheMaterialEffects
 
MuonServiceProxytheService
 
std::string theSimModuleLabel_
 
std::string theSimModuleProcess_
 
std::string theTrkModuleLabel_
 

Additional Inherited Members

- Public Types inherited from edm::EDProducer
typedef EDProducer ModuleType
 
typedef WorkerT< EDProducerWorkerType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Static Public Member Functions inherited from edm::EDProducer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDProducer
CurrentProcessingContext const * currentContext () const
 
- Protected Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Description: <one line="" class="" summary>="">

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

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

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

Definition at line 57 of file MuonSimHitProducer.h.

Constructor & Destructor Documentation

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

Definition at line 80 of file MuonSimHitProducer.cc.

References edm::hlt::Exception, edm::ParameterSet::getParameter(), edm::Service< T >::isAvailable(), MuonServiceProxy_cff::MuonServiceProxy, random, readParameters(), and theService.

80  :
81  theEstimator(iConfig.getParameter<double>("Chi2EstimatorCut")),
83  {
84 
85  //
86  // Initialize the random number generator service
87  //
89  if ( ! rng.isAvailable() ) {
90  throw cms::Exception("Configuration") <<
91  "MuonSimHitProducer requires the RandomGeneratorService \n"
92  "which is not present in the configuration file. \n"
93  "You must add the service in the configuration file\n"
94  "or remove the module that requires it.";
95  }
96 
97  random = new RandomEngine(&(*rng));
98 
99  // Read relevant parameters
101  iConfig.getParameter<edm::ParameterSet>("TRACKS"),
102  iConfig.getParameter<edm::ParameterSet>("MaterialEffectsForMuons"));
103 
104  //
105  // register your products ... need to declare at least one possible product...
106  //
107  produces<edm::PSimHitContainer>("MuonCSCHits");
108  produces<edm::PSimHitContainer>("MuonDTHits");
109  produces<edm::PSimHitContainer>("MuonRPCHits");
110 
111  edm::ParameterSet serviceParameters =
112  iConfig.getParameter<edm::ParameterSet>("ServiceParameters");
113  theService = new MuonServiceProxy(serviceParameters);
114 }
T getParameter(std::string const &) const
Propagator * propagatorWithoutMaterial
bool isAvailable() const
Definition: Service.h:47
const RandomEngine * random
MuonServiceProxy * theService
void readParameters(const edm::ParameterSet &, const edm::ParameterSet &, const edm::ParameterSet &)
Chi2MeasurementEstimator theEstimator
MuonSimHitProducer::~MuonSimHitProducer ( )

Definition at line 147 of file MuonSimHitProducer.cc.

References propagatorWithoutMaterial, random, and theMaterialEffects.

148 {
149  // do anything here that needs to be done at destruction time
150  // (e.g. close files, deallocate resources etc.)
151 
152  if ( random ) {
153  delete random;
154  }
155 
158 }
MaterialEffects * theMaterialEffects
Propagator * propagatorWithoutMaterial
const RandomEngine * random

Member Function Documentation

void MuonSimHitProducer::applyMaterialEffects ( TrajectoryStateOnSurface tsosWithdEdx,
TrajectoryStateOnSurface tsos,
double  radPath 
)
private

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

Definition at line 580 of file MuonSimHitProducer.cc.

References DeDxDiscriminatorTools::charge(), TrajectoryStateOnSurface::charge(), EnergyLossSimulator::deltaMom(), MaterialEffects::energyLossSimulator(), TrajectoryStateOnSurface::globalMomentum(), TrajectoryStateOnSurface::globalPosition(), PV3DBase< T, PVType, FrameType >::mag2(), magfield, RawParticle::momentum(), EnergyLossSimulator::mostLikelyLoss(), RPCpg::mu, MaterialEffects::multipleScatteringSimulator(), MaterialEffects::muonBremsstrahlungSimulator(), position, RawParticle::setID(), MaterialEffectsSimulator::setNormalVector(), RawParticle::setVertex(), mathSSE::sqrt(), TrajectoryStateOnSurface::surface(), Surface::tangentPlane(), theMaterialEffects, MaterialEffectsSimulator::updateState(), 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().

582  {
583 
584  // The energy loss simulator
586 
587  // The multiple scattering simulator
589 
590  // The Muon Bremsstrahlung simulator
592 
593 
594  // Initialize the Particle position, momentum and energy
595  const Surface& nextSurface = tsos.surface();
596  GlobalPoint gPos = energyLoss ? tsos.globalPosition() : tsosWithdEdx.globalPosition();
597  GlobalVector gMom = energyLoss ? tsos.globalMomentum() : tsosWithdEdx.globalMomentum();
598  double mu = 0.1056583692;
599  double en = std::sqrt(gMom.mag2()+mu*mu);
600 
601  // And now create the Particle
602  XYZTLorentzVector position(gPos.x(),gPos.y(),gPos.z(),0.);
603  XYZTLorentzVector momentum(gMom.x(),gMom.y(),gMom.z(),en);
604  float charge = (float)(tsos.charge());
605  ParticlePropagator theMuon(momentum,position,charge,0);
606  theMuon.setID(-(int)charge*13);
607 
608  // Recompute the energy loss to get the fluctuations
609  if ( energyLoss ) {
610  // Difference between with and without dE/dx (average only)
611  // (for corrections once fluctuations are applied)
612  GlobalPoint gPosWithdEdx = tsosWithdEdx.globalPosition();
613  GlobalVector gMomWithdEdx = tsosWithdEdx.globalMomentum();
614  double enWithdEdx = std::sqrt(gMomWithdEdx.mag2()+mu*mu);
616  deltaPos(gPosWithdEdx.x()-gPos.x(),gPosWithdEdx.y()-gPos.y(),
617  gPosWithdEdx.z()-gPos.z(),0.);
619  deltaMom(gMomWithdEdx.x()-gMom.x(),gMomWithdEdx.y()-gMom.y(),
620  gMomWithdEdx.z()-gMom.z(), enWithdEdx -en);
621 
622  // Energy loss in iron (+ fluctuations)
623  energyLoss->updateState(theMuon,radPath);
624 
625  // Correcting factors to account for slight differences in material descriptions
626  // (Material description is more accurate in the stepping helix propagator)
627  radPath *= -deltaMom.E()/energyLoss->mostLikelyLoss();
628  double fac = energyLoss->deltaMom().E()/energyLoss->mostLikelyLoss();
629 
630  // Particle momentum & position after energy loss + fluctuation
631  XYZTLorentzVector theNewMomentum = theMuon.momentum() + energyLoss->deltaMom() + fac * deltaMom;
632  XYZTLorentzVector theNewPosition = theMuon.vertex() + fac * deltaPos;
633  fac = (theNewMomentum.E()*theNewMomentum.E()-mu*mu)/theNewMomentum.Vect().Mag2();
634  fac = fac>0.? std::sqrt(fac) : 1E-9;
635  theMuon.SetXYZT(theNewMomentum.Px()*fac,theNewMomentum.Py()*fac,
636  theNewMomentum.Pz()*fac,theNewMomentum.E());
637  theMuon.setVertex(theNewPosition);
638 
639  }
640 
641  // Does the actual mutliple scattering
642  if ( multipleScattering ) {
643  // Pass the vector normal to the "next" surface
644  GlobalVector normal = nextSurface.tangentPlane(tsos.globalPosition())->normalVector();
645  multipleScattering->setNormalVector(normal);
646  // Compute the amount of multiple scattering after a given path length
647  multipleScattering->updateState(theMuon,radPath);
648  }
649 
650  // Muon Bremsstrahlung
651  if ( bremsstrahlung ) {
652  // Compute the amount of Muon Bremsstrahlung after given path length
653  bremsstrahlung->updateState(theMuon,radPath);
654  }
655 
656 
657  // Fill the propagated state
658  GlobalPoint propagatedPosition(theMuon.X(),theMuon.Y(),theMuon.Z());
659  GlobalVector propagatedMomentum(theMuon.Px(),theMuon.Py(),theMuon.Pz());
660  GlobalTrajectoryParameters propagatedGtp(propagatedPosition,propagatedMomentum,(int)charge,magfield);
661  tsosWithdEdx = TrajectoryStateOnSurface(propagatedGtp,nextSurface);
662 
663 }
double mostLikelyLoss() const
Return most probable energy loss.
const MagneticField * magfield
T mag2() const
Definition: PV3DBase.h:66
void setNormalVector(const GlobalVector &normal)
Sets the vector normal to the surface traversed.
MaterialEffects * theMaterialEffects
MultipleScatteringSimulator * multipleScatteringSimulator() const
Return the Multiple Scattering engine.
T y() const
Definition: PV3DBase.h:63
GlobalPoint globalPosition() const
double charge(const std::vector< uint8_t > &Ampls)
static int position[TOTALCHAMBERS][3]
Definition: ReadPGInfo.cc:509
MuonBremsstrahlungSimulator * muonBremsstrahlungSimulator() const
Return the Muon Bremsstrahlung engine.
const SurfaceType & surface() const
T sqrt(T t)
Definition: SSEVec.h:48
T z() const
Definition: PV3DBase.h:64
const XYZTLorentzVector & deltaMom() const
Returns the actual energy lost.
const int mu
Definition: Constants.h:23
virtual ReferenceCountingPointer< TangentPlane > tangentPlane(const GlobalPoint &) const =0
void updateState(ParticlePropagator &myTrack, double radlen)
Compute the material effect (calls the sub class)
GlobalVector globalMomentum() const
EnergyLossSimulator * energyLossSimulator() const
Return the Energy Loss engine.
T x() const
Definition: PV3DBase.h:62
math::XYZTLorentzVector XYZTLorentzVector
Definition: RawParticle.h:15
void MuonSimHitProducer::beginRun ( edm::Run const &  run,
const edm::EventSetup es 
)
overrideprivatevirtual

Reimplemented from edm::EDProducer.

Definition at line 118 of file MuonSimHitProducer.cc.

References Propagator::clone(), cscGeom, dtGeom, edm::EventSetup::get(), magfield, MuonServiceProxy::propagator(), LargeD0_PixelPairStep_cff::propagator, propagatorWithMaterial, propagatorWithoutMaterial, rpcGeom, theService, and MuonServiceProxy::update().

118  {
119 
120  //services
122  edm::ESHandle<DTGeometry> dtGeometry;
123  edm::ESHandle<CSCGeometry> cscGeometry;
124  edm::ESHandle<RPCGeometry> rpcGeometry;
126 
127  es.get<IdealMagneticFieldRecord>().get(magField);
128  es.get<MuonGeometryRecord>().get("MisAligned",dtGeometry);
129  es.get<MuonGeometryRecord>().get("MisAligned",cscGeometry);
130  es.get<MuonGeometryRecord>().get(rpcGeometry);
131 
132  magfield = &(*magField);
133  dtGeom = &(*dtGeometry);
134  cscGeom = &(*cscGeometry);
135  rpcGeom = &(*rpcGeometry);
136 
137  theService->update(es);
138 
139  // A few propagators
140  propagatorWithMaterial = &(*(theService->propagator("SteppingHelixPropagatorAny")));
142  SteppingHelixPropagator* SHpropagator = dynamic_cast<SteppingHelixPropagator*> (propagatorWithoutMaterial); // Beuark!
143  SHpropagator->setMaterialMode(true); // switches OFF material effects;
144 
145 }
void update(const edm::EventSetup &setup)
update the services each event
const MagneticField * magfield
const DTGeometry * dtGeom
virtual Propagator * clone() const =0
const RPCGeometry * rpcGeom
const Propagator * propagatorWithMaterial
Propagator * propagatorWithoutMaterial
const CSCGeometry * cscGeom
const T & get() const
Definition: EventSetup.h:55
MuonServiceProxy * theService
edm::ESHandle< Propagator > propagator(std::string propagatorName) const
get the propagator
void MuonSimHitProducer::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprivatevirtual

Implements edm::EDProducer.

Definition at line 168 of file MuonSimHitProducer.cc.

References abs, alongMomentum, anyDirection, applyMaterialEffects(), PV3DBase< T, PVType, FrameType >::basicVector(), Surface::bounds(), RPCGeometry::chamber(), CSCGeometry::chamber(), DTGeometry::chamber(), DTTopology::channel(), CoreSimTrack::charge(), gather_cfg::cout, Vector3DBase< T, FrameTag >::cross(), GeomDetEnumerators::CSC, cscGeom, MuonServiceProxy::detLayerGeometry(), GeomDetEnumerators::DT, dtGeom, MuonPatternRecoDumper::dumpLayer(), CSCDetId::endcap(), cmsRelvalreport::exit, create_public_lumi_plots::exp, fCSC, fDT, DTTopology::firstChannel(), RandomEngine::flatShoot(), TrajectoryStateOnSurface::freeTrajectoryState(), GeomDet::geographicalId(), edm::Event::getByLabel(), SteppingHelixStateInfo::getStateOnSurface(), TrajectoryStateOnSurface::globalMomentum(), TrajectoryStateOnSurface::globalPosition(), i, RPCGeometry::idToDetUnit(), DTGeometry::idToDetUnit(), CSCGeometry::idToDetUnit(), Bounds::inside(), CSCLayerGeometry::inside(), TrajectoryStateOnSurface::isValid(), kCSC, kDT, DTTopology::lastChannel(), DTLayerId::layer(), CSCDetId::layer(), RPCDetId::layer(), CSCChamber::layers(), TrajectoryStateOnSurface::localMomentum(), create_public_lumi_plots::log, m, PV3DBase< T, PVType, FrameType >::mag(), magfield, CoreSimTrack::momentum(), n, getHLTPrescaleColumns::path, SteppingHelixStateInfo::path(), PV3DBase< T, PVType, FrameType >::phi(), pi, evf::utils::pid, PlaneBuilder::plane(), position, Propagator::propagateWithPath(), propagatorWithMaterial, propagatorWithoutMaterial, edm::Event::put(), SteppingHelixStateInfo::radPath(), random, DetId::rawId(), relativeConstraints::ring, RPCDetId::ring(), CSCDetId::ring(), RPCDetId::roll(), RPCChamber::rolls(), makeMuonMisalignmentScenario::rot, GeomDetEnumerators::RPCBarrel, GeomDetEnumerators::RPCEndcap, rpcGeom, DTChamberId::sector(), RPCDetId::sector(), mathSSE::sqrt(), DTChamberId::station(), relativeConstraints::station, CSCDetId::station(), RPCDetId::station(), GeomDet::subDetector(), RPCDetId::subsector(), DTSuperLayerId::superlayer(), DTChamber::superLayers(), GeomDet::surface(), theEstimator, theMaterialEffects, theService, theSimModuleLabel_, theSimModuleProcess_, PV3DBase< T, PVType, FrameType >::theta(), Bounds::thickness(), GeomDet::toLocal(), SimTrack::trackerSurfaceMomentum(), SimTrack::trackerSurfacePosition(), CoreSimTrack::trackId(), TrajectoryStateOnSurface::transverseCurvature(), CoreSimTrack::type(), Vector3DBase< T, FrameTag >::unit(), SimTrack::vertIndex(), DTChamberId::wheel(), x, PV3DBase< T, PVType, FrameType >::x(), detailsBasic3DVector::y, PV3DBase< T, PVType, FrameType >::y(), detailsBasic3DVector::z, and PV3DBase< T, PVType, FrameType >::z().

Referenced by JSONExport.JsonExport::export(), HTMLExport.HTMLExport::export(), and HTMLExport.HTMLExportStatic::export().

168  {
169  // using namespace edm;
170  // using namespace std;
171 
172  MuonPatternRecoDumper dumper;
173 
176  std::vector<PSimHit> theCSCHits;
177  std::vector<PSimHit> theDTHits;
178  std::vector<PSimHit> theRPCHits;
179 
182  iEvent.getByLabel(theSimModuleLabel_,simVertices);
183 
184  for ( unsigned int itrk=0; itrk<simMuons->size(); itrk++ ) {
185  const SimTrack &mySimTrack = (*simMuons)[itrk];
186  math::XYZTLorentzVector mySimP4(mySimTrack.momentum().x(),
187  mySimTrack.momentum().y(),
188  mySimTrack.momentum().z(),
189  mySimTrack.momentum().t());
190 
191  // Decaying hadrons are now in the list, and so are their muon daughter
192  // Ignore the hadrons here.
193  int pid = mySimTrack.type();
194  if ( abs(pid) != 13 ) continue;
195 
196  double t0 = 0;
197  GlobalPoint initialPosition;
198  int ivert = mySimTrack.vertIndex();
199  if ( ivert >= 0 ) {
200  t0 = (*simVertices)[ivert].position().t();
201  GlobalPoint xyzzy((*simVertices)[ivert].position().x(),
202  (*simVertices)[ivert].position().y(),
203  (*simVertices)[ivert].position().z());
204  initialPosition = xyzzy;
205  }
206 //
207 // Presumably t0 has dimensions of cm if not mm?
208 // Convert to ns for internal calculations.
209 // I wonder where we should get c from?
210 //
211  double tof = t0/29.98;
212 
213 #ifdef FAMOS_DEBUG
214  std::cout << " ===> MuonSimHitProducer::reconstruct() found SIMTRACK - pid = "
215  << pid ;
216  std::cout << " : pT = " << mySimP4.Pt()
217  << ", eta = " << mySimP4.Eta()
218  << ", phi = " << mySimP4.Phi() << std::endl;
219 #endif
220 
221 //
222 // Produce muons sim hits starting from undecayed simulated muons
223 //
224 
225  GlobalPoint startingPosition(mySimTrack.trackerSurfacePosition().x(),
226  mySimTrack.trackerSurfacePosition().y(),
227  mySimTrack.trackerSurfacePosition().z());
228  GlobalVector startingMomentum(mySimTrack.trackerSurfaceMomentum().x(),
229  mySimTrack.trackerSurfaceMomentum().y(),
230  mySimTrack.trackerSurfaceMomentum().z());
231 //
232 // Crap... there's no time-of-flight to the trackerSurfacePosition()...
233 // So, this will be wrong when the curvature can't be neglected, but that
234 // will be rather seldom... May as well ignore the mass too.
235 //
236  GlobalVector dtracker = startingPosition-initialPosition;
237  tof += dtracker.mag()/29.98;
238 
239 #ifdef FAMOS_DEBUG
240  std::cout << " the Muon START position " << startingPosition << std::endl;
241  std::cout << " the Muon START momentum " << startingMomentum << std::endl;
242 #endif
243 
244 //
245 // Some magic to define a TrajectoryStateOnSurface
246 //
247  PlaneBuilder pb;
248  GlobalVector zAxis = startingMomentum.unit();
249  GlobalVector yAxis(zAxis.y(),-zAxis.x(),0);
250  GlobalVector xAxis = yAxis.cross(zAxis);
251  Surface::RotationType rot = Surface::RotationType(xAxis,yAxis,zAxis);
252  PlaneBuilder::ReturnType startingPlane = pb.plane(startingPosition,rot);
253  GlobalTrajectoryParameters gtp(startingPosition,
254  startingMomentum,
255  (int)mySimTrack.charge(),
256  magfield);
257  TrajectoryStateOnSurface startingState(gtp,*startingPlane);
258 
259  std::vector<const DetLayer *> navLayers;
260  if ( fabs(startingState.globalMomentum().eta()) > 4.5 ) {
261  navLayers = navigation.compatibleEndcapLayers(*(startingState.freeState()),
262  alongMomentum);
263  }
264  else {
265  navLayers = navigation.compatibleLayers(*(startingState.freeState()),
266  alongMomentum);
267  }
268  /*
269  edm::ESHandle<Propagator> propagator =
270  theService->propagator("SteppingHelixPropagatorAny");
271  */
272 
273  if ( navLayers.empty() ) continue;
274 
275 #ifdef FAMOS_DEBUG
276  std::cout << "Found " << navLayers.size()
277  << " compatible DetLayers..." << std::endl;
278 #endif
279 
280  TrajectoryStateOnSurface propagatedState = startingState;
281  for ( unsigned int ilayer=0; ilayer<navLayers.size(); ilayer++ ) {
282 
283 #ifdef FAMOS_DEBUG
284  std::cout << "Propagating to layer " << ilayer << " "
285  << dumper.dumpLayer(navLayers[ilayer])
286  << std::endl;
287 #endif
288 
289  std::vector<DetWithState> comps =
290  navLayers[ilayer]->compatibleDets(propagatedState,*propagatorWithMaterial,theEstimator);
291  if ( comps.empty() ) continue;
292 
293 #ifdef FAMOS_DEBUG
294  std::cout << "Propagating " << propagatedState << std::endl;
295 #endif
296 
297  // Starting momentum
298  double pi = propagatedState.globalMomentum().mag();
299 
300  // Propagate with material effects (dE/dx average only)
301  SteppingHelixStateInfo shsStart(*(propagatedState.freeTrajectoryState()));
302  const SteppingHelixStateInfo& shsDest =
303  ((const SteppingHelixPropagator*)propagatorWithMaterial)->propagate(shsStart,navLayers[ilayer]->surface());
304  std::pair<TrajectoryStateOnSurface,double> next(shsDest.getStateOnSurface(navLayers[ilayer]->surface()),
305  shsDest.path());
306  // No need to continue if there is no valid propagation available.
307  // This happens rarely (~0.1% of ttbar events)
308  if ( !next.first.isValid() ) continue;
309  // This is the estimate of the number of radiation lengths traversed,
310  // together with the total path length
311  double radPath = shsDest.radPath();
312  double pathLength = next.second;
313 
314  // Now propagate without dE/dx (average)
315  // [To add the dE/dx fluctuations to the actual dE/dx]
316  std::pair<TrajectoryStateOnSurface,double> nextNoMaterial =
317  propagatorWithoutMaterial->propagateWithPath(propagatedState,navLayers[ilayer]->surface());
318 
319  // Update the propagated state
320  propagatedState = next.first;
321  double pf = propagatedState.globalMomentum().mag();
322 
323  // Insert dE/dx fluctuations and multiple scattering
324  // Skip this step if nextNoMaterial.first is not valid
325  // This happens rarely (~0.02% of ttbar events)
326  if ( theMaterialEffects && nextNoMaterial.first.isValid() ) applyMaterialEffects(propagatedState,nextNoMaterial.first,radPath);
327  // Check that the 'shaken' propagatedState is still valid, otherwise continue
328  if ( !propagatedState.isValid() ) continue;
329  // (No evidence that this ever happens)
330 //
331 // Consider this... 1 GeV muon has a velocity that is only 0.5% slower than c...
332 // We probably can safely ignore the mass for anything that makes it out to the
333 // muon chambers.
334 //
335  double pavg = 0.5*(pi+pf);
336  double m = mySimP4.M();
337  double rbeta = sqrt(1+m*m/(pavg*pavg))/29.98;
338  double dtof = pathLength*rbeta;
339 
340 #ifdef FAMOS_DEBUG
341  std::cout << "Propagated to next surface... path length = " << pathLength
342  << " cm, dTOF = " << dtof << " ns" << std::endl;
343 #endif
344 
345  tof += dtof;
346 
347  for ( unsigned int icomp=0; icomp<comps.size(); icomp++ )
348  {
349  const GeomDet *gd = comps[icomp].first;
350  if ( gd->subDetector() == GeomDetEnumerators::DT ) {
351  DTChamberId id(gd->geographicalId());
352  const DTChamber *chamber = dtGeom->chamber(id);
353  std::vector<const DTSuperLayer *> superlayer = chamber->superLayers();
354  for ( unsigned int isl=0; isl<superlayer.size(); isl++ ) {
355  std::vector<const DTLayer *> layer = superlayer[isl]->layers();
356  for ( unsigned int ilayer=0; ilayer<layer.size(); ilayer++ ) {
357  DTLayerId lid = layer[ilayer]->id();
358 #ifdef FAMOS_DEBUG
359  std::cout << " Extrapolated to DT ("
360  << lid.wheel() << ","
361  << lid.station() << ","
362  << lid.sector() << ","
363  << lid.superlayer() << ","
364  << lid.layer() << ")" << std::endl;
365 #endif
366 
367  const GeomDetUnit *det = dtGeom->idToDetUnit(lid);
368 
369  HelixArbitraryPlaneCrossing crossing(propagatedState.globalPosition().basicVector(),
370  propagatedState.globalMomentum().basicVector(),
371  propagatedState.transverseCurvature(),
372  anyDirection);
373  std::pair<bool,double> path = crossing.pathLength(det->surface());
374  if ( ! path.first ) continue;
375  LocalPoint lpos = det->toLocal(GlobalPoint(crossing.position(path.second)));
376  if ( ! det->surface().bounds().inside(lpos) ) continue;
377  const DTTopology& dtTopo = layer[ilayer]->specificTopology();
378  int wire = dtTopo.channel(lpos);
379  if (wire - dtTopo.firstChannel() < 0 || wire - dtTopo.lastChannel() > 0 ) continue;
380  // no drift cell here (on the chamber edge or just outside)
381  // this hit would otherwise be discarded downstream in the digitizer
382 
383  DTWireId wid(lid,wire);
384  double thickness = det->surface().bounds().thickness();
385  LocalVector lmom = det->toLocal(GlobalVector(crossing.direction(path.second)));
386  lmom = lmom.unit()*propagatedState.localMomentum().mag();
387 
388  // Factor that takes into account the (rec)hits lost because of delta's, etc.:
389  // (Not fully satisfactory patch, but it seems to work...)
390  double pmu = lmom.mag();
391  double theDTHitIneff = pmu>0? exp(kDT*log(pmu)+fDT):0.;
392  if (random->flatShoot()<theDTHitIneff) continue;
393 
394  double eloss = 0;
395  double pz = fabs(lmom.z());
396  LocalPoint entry = lpos - 0.5*thickness*lmom/pz;
397  LocalPoint exit = lpos + 0.5*thickness*lmom/pz;
398  double dtof = path.second*rbeta;
399  int trkid = mySimTrack.trackId();
400  unsigned int id = wid.rawId();
401  short unsigned int processType = 2;
402  PSimHit hit(entry,exit,lmom.mag(),
403  tof+dtof,eloss,pid,id,trkid,lmom.theta(),lmom.phi(),processType);
404  theDTHits.push_back(hit);
405 
406  }
407  }
408  }
409  else if ( gd->subDetector() == GeomDetEnumerators::CSC ) {
410  CSCDetId id(gd->geographicalId());
411  const CSCChamber *chamber = cscGeom->chamber(id);
412  std::vector<const CSCLayer *> layers = chamber->layers();
413  for ( unsigned int ilayer=0; ilayer<layers.size(); ilayer++ ) {
414  CSCDetId lid = layers[ilayer]->id();
415 
416 #ifdef FAMOS_DEBUG
417  std::cout << " Extrapolated to CSC ("
418  << lid.endcap() << ","
419  << lid.ring() << ","
420  << lid.station() << ","
421  << lid.layer() << ")" << std::endl;
422 #endif
423 
424  const GeomDetUnit *det = cscGeom->idToDetUnit(lid);
425  HelixArbitraryPlaneCrossing crossing(propagatedState.globalPosition().basicVector(),
426  propagatedState.globalMomentum().basicVector(),
427  propagatedState.transverseCurvature(),
428  anyDirection);
429  std::pair<bool,double> path = crossing.pathLength(det->surface());
430  if ( ! path.first ) continue;
431  LocalPoint lpos = det->toLocal(GlobalPoint(crossing.position(path.second)));
432  // For CSCs the Bounds are for chamber frames not gas regions
433  // if ( ! det->surface().bounds().inside(lpos) ) continue;
434  // New function knows where the 'active' volume is:
435  const CSCLayerGeometry* laygeom = layers[ilayer]->geometry();
436  if ( ! laygeom->inside( lpos ) ) continue;
437  //double thickness = laygeom->thickness(); gives number which is about 20 times too big
438  double thickness = det->surface().bounds().thickness(); // this one works much better...
439  LocalVector lmom = det->toLocal(GlobalVector(crossing.direction(path.second)));
440  lmom = lmom.unit()*propagatedState.localMomentum().mag();
441 
442  // Factor that takes into account the (rec)hits lost because of delta's, etc.:
443  // (Not fully satisfactory patch, but it seems to work...)
444  double pmu = lmom.mag();
445  double theCSCHitIneff = pmu>0? exp(kCSC*log(pmu)+fCSC):0.;
446  // Take into account the different geometry in ME11:
447  if (id.station()==1 && id.ring()==1) theCSCHitIneff = theCSCHitIneff*0.442;
448  if (random->flatShoot()<theCSCHitIneff) continue;
449 
450  double eloss = 0;
451  double pz = fabs(lmom.z());
452  LocalPoint entry = lpos - 0.5*thickness*lmom/pz;
453  LocalPoint exit = lpos + 0.5*thickness*lmom/pz;
454  double dtof = path.second*rbeta;
455  int trkid = mySimTrack.trackId();
456  unsigned int id = lid.rawId();
457  short unsigned int processType = 2;
458  PSimHit hit(entry,exit,lmom.mag(),
459  tof+dtof,eloss,pid,id,trkid,lmom.theta(),lmom.phi(),processType);
460  theCSCHits.push_back(hit);
461  }
462  }
463  else if ( gd->subDetector() == GeomDetEnumerators::RPCBarrel ||
465  RPCDetId id(gd->geographicalId());
466  const RPCChamber *chamber = rpcGeom->chamber(id);
467  std::vector<const RPCRoll *> roll = chamber->rolls();
468  for ( unsigned int iroll=0; iroll<roll.size(); iroll++ ) {
469  RPCDetId rid = roll[iroll]->id();
470 
471 #ifdef FAMOS_DEBUG
472  std::cout << " Extrapolated to RPC ("
473  << rid.ring() << ","
474  << rid.station() << ","
475  << rid.sector() << ","
476  << rid.subsector() << ","
477  << rid.layer() << ","
478  << rid.roll() << ")" << std::endl;
479 #endif
480 
481  const GeomDetUnit *det = rpcGeom->idToDetUnit(rid);
482  HelixArbitraryPlaneCrossing crossing(propagatedState.globalPosition().basicVector(),
483  propagatedState.globalMomentum().basicVector(),
484  propagatedState.transverseCurvature(),
485  anyDirection);
486  std::pair<bool,double> path = crossing.pathLength(det->surface());
487  if ( ! path.first ) continue;
488  LocalPoint lpos = det->toLocal(GlobalPoint(crossing.position(path.second)));
489  if ( ! det->surface().bounds().inside(lpos) ) continue;
490  double thickness = det->surface().bounds().thickness();
491  LocalVector lmom = det->toLocal(GlobalVector(crossing.direction(path.second)));
492  lmom = lmom.unit()*propagatedState.localMomentum().mag();
493  double eloss = 0;
494  double pz = fabs(lmom.z());
495  LocalPoint entry = lpos - 0.5*thickness*lmom/pz;
496  LocalPoint exit = lpos + 0.5*thickness*lmom/pz;
497  double dtof = path.second*rbeta;
498  int trkid = mySimTrack.trackId();
499  unsigned int id = rid.rawId();
500  short unsigned int processType = 2;
501  PSimHit hit(entry,exit,lmom.mag(),
502  tof+dtof,eloss,pid,id,trkid,lmom.theta(),lmom.phi(),processType);
503  theRPCHits.push_back(hit);
504  }
505  }
506  else {
507  std::cout << "Extrapolated to unknown subdetector '" << gd->subDetector() << "'..." << std::endl;
508  }
509  }
510  }
511  }
512 
513  std::auto_ptr<edm::PSimHitContainer> pcsc(new edm::PSimHitContainer);
514  int n = 0;
515  for ( std::vector<PSimHit>::const_iterator i = theCSCHits.begin();
516  i != theCSCHits.end(); i++ ) {
517  pcsc->push_back(*i);
518  n += 1;
519  }
520  iEvent.put(pcsc,"MuonCSCHits");
521 
522  std::auto_ptr<edm::PSimHitContainer> pdt(new edm::PSimHitContainer);
523  n = 0;
524  for ( std::vector<PSimHit>::const_iterator i = theDTHits.begin();
525  i != theDTHits.end(); i++ ) {
526  pdt->push_back(*i);
527  n += 1;
528  }
529  iEvent.put(pdt,"MuonDTHits");
530 
531  std::auto_ptr<edm::PSimHitContainer> prpc(new edm::PSimHitContainer);
532  n = 0;
533  for ( std::vector<PSimHit>::const_iterator i = theRPCHits.begin();
534  i != theRPCHits.end(); i++ ) {
535  prpc->push_back(*i);
536  n += 1;
537  }
538  iEvent.put(prpc,"MuonRPCHits");
539 
540 }
edm::ESHandle< MuonDetLayerGeometry > detLayerGeometry() const
get the detLayer geometry
const math::XYZVectorD & trackerSurfacePosition() const
Definition: SimTrack.h:36
const MagneticField * magfield
int i
Definition: DBlmapReader.cc:9
const std::vector< const CSCLayer * > & layers() const
Return all layers.
Definition: CSCChamber.h:57
const DTGeometry * dtGeom
std::string dumpLayer(const DetLayer *layer) const
virtual bool inside(const Local3DPoint &) const =0
Determine if the point is inside the bounds.
virtual const GeomDetUnit * idToDetUnit(DetId) const
Return the pointer to the GeomDetUnit corresponding to a given DetId.
Definition: RPCGeometry.cc:48
FreeTrajectoryState * freeTrajectoryState(bool withErrors=true) const
int channel(const LocalPoint &p) const
Definition: DTTopology.cc:81
MaterialEffects * theMaterialEffects
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
ReturnType plane(const PositionType &pos, const RotationType &rot) const
Definition: PlaneBuilder.h:22
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
T y() const
Definition: PV3DBase.h:63
virtual SubDetector subDetector() const =0
Which subdetector.
#define abs(x)
Definition: mlp_lapack.h:159
GlobalPoint globalPosition() const
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
Definition: GeomDet.h:62
const Bounds & bounds() const
Definition: Surface.h:128
int layer() const
Return the layer number.
Definition: DTLayerId.h:55
float charge() const
charge
Definition: CoreSimTrack.cc:3
int layer() const
Definition: CSCDetId.h:63
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:35
float float float z
static int position[TOTALCHAMBERS][3]
Definition: ReadPGInfo.cc:509
int firstChannel() const
Returns the wire number of the first wire.
Definition: DTTopology.h:80
Geom::Theta< T > theta() const
Definition: PV3DBase.h:75
const RPCGeometry * rpcGeom
uint32_t rawId() const
get the raw id
Definition: DetId.h:45
LocalVector localMomentum() const
int endcap() const
Definition: CSCDetId.h:95
virtual float thickness() const =0
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:30
int lastChannel() const
Returns the wire number of the last wire.
Definition: DTTopology.h:82
const RPCChamber * chamber(RPCDetId id) const
Definition: RPCGeometry.cc:71
virtual std::pair< TrajectoryStateOnSurface, double > propagateWithPath(const FreeTrajectoryState &, const Surface &) const
Definition: Propagator.cc:77
T mag() const
Definition: PV3DBase.h:67
int roll() const
Definition: RPCDetId.h:123
std::string theSimModuleProcess_
virtual const GeomDetUnit * idToDetUnit(DetId) const
Return the pointer to the GeomDetUnit corresponding to a given DetId.
Definition: CSCGeometry.cc:93
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:94
const std::vector< const DTSuperLayer * > & superLayers() const
Return the superlayers in the chamber.
Definition: DTChamber.cc:62
int ring() const
Definition: RPCDetId.h:75
const Propagator * propagatorWithMaterial
TrajectoryStateOnSurface getStateOnSurface(const Surface &surf, bool returnTangentPlane=false) const
T sqrt(T t)
Definition: SSEVec.h:48
Propagator * propagatorWithoutMaterial
std::pair< std::string, MonitorElement * > entry
Definition: ME_MAP.h:8
Vector3DBase< typename PreciseFloatType< T, U >::Type, FrameTag > cross(const Vector3DBase< U, FrameTag > &v) const
Definition: Vector3DBase.h:119
T z() const
Definition: PV3DBase.h:64
const std::vector< const RPCRoll * > & rolls() const
Return the Rolls.
Definition: RPCChamber.cc:70
const CSCGeometry * cscGeom
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:72
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
int vertIndex() const
index of the vertex in the Event container (-1 if no vertex)
Definition: SimTrack.h:29
int superlayer() const
Return the superlayer number (deprecated method name)
int ring() const
Definition: CSCDetId.h:77
Vector3DBase unit() const
Definition: Vector3DBase.h:57
int layer() const
Definition: RPCDetId.h:111
unsigned int trackId() const
Definition: CoreSimTrack.h:49
const math::XYZTLorentzVectorD & trackerSurfaceMomentum() const
Definition: SimTrack.h:38
const CSCChamber * chamber(CSCDetId id) const
Return the chamber corresponding to given DetId.
Definition: CSCGeometry.cc:118
const RandomEngine * random
const DTChamber * chamber(DTChamberId id) const
Return a DTChamber given its id.
Definition: DTGeometry.cc:102
double flatShoot(double xmin=0.0, double xmax=1.0) const
Definition: RandomEngine.h:30
int sector() const
Sector id: the group of chambers at same phi (and increasing r)
Definition: RPCDetId.h:105
MuonServiceProxy * theService
void applyMaterialEffects(TrajectoryStateOnSurface &tsosWithdEdx, TrajectoryStateOnSurface &tsos, double radPath)
Simulate material effects in iron (dE/dx, multiple scattering)
int subsector() const
SubSector id : some sectors are divided along the phi direction in subsectors (from 1 to 4 in Barrel...
Definition: RPCDetId.h:117
int type() const
particle type (HEP PDT convension)
Definition: CoreSimTrack.h:40
std::string theSimModuleLabel_
GlobalVector globalMomentum() const
int sector() const
Definition: DTChamberId.h:63
const math::XYZTLorentzVectorD & momentum() const
particle info...
Definition: CoreSimTrack.h:36
TkRotation< float > RotationType
virtual const GeomDetUnit * idToDetUnit(DetId) const
Return the pointer to the GeomDetUnit corresponding to a given DetId.
Definition: DTGeometry.cc:73
bool inside(const Local3DPoint &, const LocalError &, float scale=1.f) const
int station() const
Definition: CSCDetId.h:88
tuple cout
Definition: gather_cfg.py:121
std::vector< PSimHit > PSimHitContainer
double pi
Definition: DDAxes.h:10
int station() const
Return the station number.
Definition: DTChamberId.h:53
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:47
T x() const
Definition: PV3DBase.h:62
const BasicVectorType & basicVector() const
Definition: PV3DBase.h:56
Global3DVector GlobalVector
Definition: GlobalVector.h:10
Chi2MeasurementEstimator theEstimator
int station() const
Definition: RPCDetId.h:99
void MuonSimHitProducer::readParameters ( const edm::ParameterSet fastMuons,
const edm::ParameterSet fastTracks,
const edm::ParameterSet matEff 
)
private

Definition at line 543 of file MuonSimHitProducer.cc.

References fCSC, fDT, fullPattern_, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), kCSC, kDT, random, AlCaHLTBitMon_QueryRunRegistry::string, theMaterialEffects, theSimModuleLabel_, theSimModuleProcess_, and theTrkModuleLabel_.

Referenced by MuonSimHitProducer().

545  {
546  // Muons
547  theSimModuleLabel_ = fastMuons.getParameter<std::string>("simModuleLabel");
548  theSimModuleProcess_ = fastMuons.getParameter<std::string>("simModuleProcess");
549  theTrkModuleLabel_ = fastMuons.getParameter<std::string>("trackModuleLabel");
550  std::vector<double> simHitIneffDT = fastMuons.getParameter<std::vector<double> >("simHitDTIneffParameters");
551  std::vector<double> simHitIneffCSC = fastMuons.getParameter<std::vector<double> >("simHitCSCIneffParameters");
552  kDT = simHitIneffDT[0];
553  fDT = simHitIneffDT[1];
554  kCSC = simHitIneffCSC[0];
555  fCSC = simHitIneffCSC[1];
556 
557  // Tracks
558  fullPattern_ = fastTracks.getUntrackedParameter<bool>("FullPatternRecognition");
559 
560  // The following should be on LogInfo
561  // std::cout << " MUON SIM HITS: FastSimulation parameters " << std::endl;
562  // std::cout << " ============================================== " << std::endl;
563  // if ( fullPattern_ )
564  // std::cout << " The FULL pattern recognition option is turned ON" << std::endl;
565  // else
566  // std::cout << " The FAST tracking option is turned ON" << std::endl;
567 
568  // Material Effects
569  theMaterialEffects = 0;
570  if ( matEff.getParameter<bool>("PairProduction") ||
571  matEff.getParameter<bool>("Bremsstrahlung") ||
572  matEff.getParameter<bool>("MuonBremsstrahlung") ||
573  matEff.getParameter<bool>("EnergyLoss") ||
574  matEff.getParameter<bool>("MultipleScattering") )
575  theMaterialEffects = new MaterialEffects(matEff,random);
576 
577 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
MaterialEffects * theMaterialEffects
std::string theSimModuleProcess_
const RandomEngine * random
std::string theSimModuleLabel_
std::string theTrkModuleLabel_

Member Data Documentation

const CSCGeometry* MuonSimHitProducer::cscGeom
private

Definition at line 71 of file MuonSimHitProducer.h.

Referenced by beginRun(), and produce().

bool MuonSimHitProducer::doGL_
private

Definition at line 98 of file MuonSimHitProducer.h.

bool MuonSimHitProducer::doL1_
private

Definition at line 98 of file MuonSimHitProducer.h.

bool MuonSimHitProducer::doL3_
private

Definition at line 98 of file MuonSimHitProducer.h.

const DTGeometry* MuonSimHitProducer::dtGeom
private

Definition at line 70 of file MuonSimHitProducer.h.

Referenced by beginRun(), and produce().

double MuonSimHitProducer::fCSC
private

Definition at line 88 of file MuonSimHitProducer.h.

Referenced by produce(), and readParameters().

double MuonSimHitProducer::fDT
private

Definition at line 86 of file MuonSimHitProducer.h.

Referenced by produce(), and readParameters().

bool MuonSimHitProducer::fullPattern_
private

Definition at line 97 of file MuonSimHitProducer.h.

Referenced by readParameters().

double MuonSimHitProducer::kCSC
private

Definition at line 87 of file MuonSimHitProducer.h.

Referenced by produce(), and readParameters().

double MuonSimHitProducer::kDT
private

Definition at line 85 of file MuonSimHitProducer.h.

Referenced by produce(), and readParameters().

const MagneticField* MuonSimHitProducer::magfield
private

Definition at line 69 of file MuonSimHitProducer.h.

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

const Propagator* MuonSimHitProducer::propagatorWithMaterial
private

Definition at line 73 of file MuonSimHitProducer.h.

Referenced by beginRun(), and produce().

Propagator* MuonSimHitProducer::propagatorWithoutMaterial
private

Definition at line 74 of file MuonSimHitProducer.h.

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

const RandomEngine* MuonSimHitProducer::random
private
const RPCGeometry* MuonSimHitProducer::rpcGeom
private

Definition at line 72 of file MuonSimHitProducer.h.

Referenced by beginRun(), and produce().

Chi2MeasurementEstimator MuonSimHitProducer::theEstimator
private

Definition at line 67 of file MuonSimHitProducer.h.

Referenced by produce().

MaterialEffects* MuonSimHitProducer::theMaterialEffects
private
MuonServiceProxy* MuonSimHitProducer::theService
private

Definition at line 66 of file MuonSimHitProducer.h.

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

std::string MuonSimHitProducer::theSimModuleLabel_
private

Definition at line 99 of file MuonSimHitProducer.h.

Referenced by produce(), and readParameters().

std::string MuonSimHitProducer::theSimModuleProcess_
private

Definition at line 99 of file MuonSimHitProducer.h.

Referenced by produce(), and readParameters().

std::string MuonSimHitProducer::theTrkModuleLabel_
private

Definition at line 99 of file MuonSimHitProducer.h.

Referenced by readParameters().