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 ()
 
ModuleDescription const & moduleDescription () const
 
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 ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Member Functions

void applyMaterialEffects (TrajectoryStateOnSurface &tsosWithdEdx, TrajectoryStateOnSurface &tsos, double radPath, RandomEngineAndDistribution const *)
 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 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
 
- 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::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 78 of file MuonSimHitProducer.cc.

References edm::ParameterSet::getParameter(), MuonServiceProxy_cff::MuonServiceProxy, readParameters(), and theService.

78  :
79  theEstimator(iConfig.getParameter<double>("Chi2EstimatorCut")),
81  {
82 
83  // Read relevant parameters
85  iConfig.getParameter<edm::ParameterSet>("TRACKS"),
86  iConfig.getParameter<edm::ParameterSet>("MaterialEffectsForMuons"));
87 
88  //
89  // register your products ... need to declare at least one possible product...
90  //
91  produces<edm::PSimHitContainer>("MuonCSCHits");
92  produces<edm::PSimHitContainer>("MuonDTHits");
93  produces<edm::PSimHitContainer>("MuonRPCHits");
94 
95  edm::ParameterSet serviceParameters =
96  iConfig.getParameter<edm::ParameterSet>("ServiceParameters");
97  theService = new MuonServiceProxy(serviceParameters);
98 }
T getParameter(std::string const &) const
Propagator * propagatorWithoutMaterial
MuonServiceProxy * theService
void readParameters(const edm::ParameterSet &, const edm::ParameterSet &, const edm::ParameterSet &)
Chi2MeasurementEstimator theEstimator
MuonSimHitProducer::~MuonSimHitProducer ( )

Definition at line 131 of file MuonSimHitProducer.cc.

References propagatorWithoutMaterial, and theMaterialEffects.

132 {
133  // do anything here that needs to be done at destruction time
134  // (e.g. close files, deallocate resources etc.)
135 
138 }
MaterialEffects * theMaterialEffects
Propagator * propagatorWithoutMaterial

Member Function Documentation

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

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

Definition at line 561 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().

564  {
565 
566  // The energy loss simulator
568 
569  // The multiple scattering simulator
571 
572  // The Muon Bremsstrahlung simulator
574 
575 
576  // Initialize the Particle position, momentum and energy
577  const Surface& nextSurface = tsos.surface();
578  GlobalPoint gPos = energyLoss ? tsos.globalPosition() : tsosWithdEdx.globalPosition();
579  GlobalVector gMom = energyLoss ? tsos.globalMomentum() : tsosWithdEdx.globalMomentum();
580  double mu = 0.1056583692;
581  double en = std::sqrt(gMom.mag2()+mu*mu);
582 
583  // And now create the Particle
584  XYZTLorentzVector position(gPos.x(),gPos.y(),gPos.z(),0.);
585  XYZTLorentzVector momentum(gMom.x(),gMom.y(),gMom.z(),en);
586  float charge = (float)(tsos.charge());
587  ParticlePropagator theMuon(momentum,position,charge,0);
588  theMuon.setID(-(int)charge*13);
589 
590  // Recompute the energy loss to get the fluctuations
591  if ( energyLoss ) {
592  // Difference between with and without dE/dx (average only)
593  // (for corrections once fluctuations are applied)
594  GlobalPoint gPosWithdEdx = tsosWithdEdx.globalPosition();
595  GlobalVector gMomWithdEdx = tsosWithdEdx.globalMomentum();
596  double enWithdEdx = std::sqrt(gMomWithdEdx.mag2()+mu*mu);
598  deltaPos(gPosWithdEdx.x()-gPos.x(),gPosWithdEdx.y()-gPos.y(),
599  gPosWithdEdx.z()-gPos.z(),0.);
601  deltaMom(gMomWithdEdx.x()-gMom.x(),gMomWithdEdx.y()-gMom.y(),
602  gMomWithdEdx.z()-gMom.z(), enWithdEdx -en);
603 
604  // Energy loss in iron (+ fluctuations)
605  energyLoss->updateState(theMuon, radPath, random);
606 
607  // Correcting factors to account for slight differences in material descriptions
608  // (Material description is more accurate in the stepping helix propagator)
609  radPath *= -deltaMom.E()/energyLoss->mostLikelyLoss();
610  double fac = energyLoss->deltaMom().E()/energyLoss->mostLikelyLoss();
611 
612  // Particle momentum & position after energy loss + fluctuation
613  XYZTLorentzVector theNewMomentum = theMuon.momentum() + energyLoss->deltaMom() + fac * deltaMom;
614  XYZTLorentzVector theNewPosition = theMuon.vertex() + fac * deltaPos;
615  fac = (theNewMomentum.E()*theNewMomentum.E()-mu*mu)/theNewMomentum.Vect().Mag2();
616  fac = fac>0.? std::sqrt(fac) : 1E-9;
617  theMuon.SetXYZT(theNewMomentum.Px()*fac,theNewMomentum.Py()*fac,
618  theNewMomentum.Pz()*fac,theNewMomentum.E());
619  theMuon.setVertex(theNewPosition);
620 
621  }
622 
623  // Does the actual mutliple scattering
624  if ( multipleScattering ) {
625  // Pass the vector normal to the "next" surface
626  GlobalVector normal = nextSurface.tangentPlane(tsos.globalPosition())->normalVector();
627  multipleScattering->setNormalVector(normal);
628  // Compute the amount of multiple scattering after a given path length
629  multipleScattering->updateState(theMuon, radPath, random);
630  }
631 
632  // Muon Bremsstrahlung
633  if ( bremsstrahlung ) {
634  // Compute the amount of Muon Bremsstrahlung after given path length
635  bremsstrahlung->updateState(theMuon, radPath, random);
636  }
637 
638 
639  // Fill the propagated state
640  GlobalPoint propagatedPosition(theMuon.X(),theMuon.Y(),theMuon.Z());
641  GlobalVector propagatedMomentum(theMuon.Px(),theMuon.Py(),theMuon.Pz());
642  GlobalTrajectoryParameters propagatedGtp(propagatedPosition,propagatedMomentum,(int)charge,magfield);
643  tsosWithdEdx = TrajectoryStateOnSurface(propagatedGtp,nextSurface);
644 
645 }
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
void updateState(ParticlePropagator &myTrack, double radlen, RandomEngineAndDistribution const *)
Compute the material effect (calls the sub class)
TRandom random
Definition: MVATrainer.cc:138
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:22
virtual ReferenceCountingPointer< TangentPlane > tangentPlane(const GlobalPoint &) const =0
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 102 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().

102  {
103 
104  //services
106  edm::ESHandle<DTGeometry> dtGeometry;
107  edm::ESHandle<CSCGeometry> cscGeometry;
108  edm::ESHandle<RPCGeometry> rpcGeometry;
110 
111  es.get<IdealMagneticFieldRecord>().get(magField);
112  es.get<MuonGeometryRecord>().get("MisAligned",dtGeometry);
113  es.get<MuonGeometryRecord>().get("MisAligned",cscGeometry);
114  es.get<MuonGeometryRecord>().get(rpcGeometry);
115 
116  magfield = &(*magField);
117  dtGeom = &(*dtGeometry);
118  cscGeom = &(*cscGeometry);
119  rpcGeom = &(*rpcGeometry);
120 
121  theService->update(es);
122 
123  // A few propagators
124  propagatorWithMaterial = &(*(theService->propagator("SteppingHelixPropagatorAny")));
126  SteppingHelixPropagator* SHpropagator = dynamic_cast<SteppingHelixPropagator*> (propagatorWithoutMaterial); // Beuark!
127  SHpropagator->setMaterialMode(true); // switches OFF material effects;
128 
129 }
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 148 of file MuonSimHitProducer.cc.

References funct::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, CSCDetId::endcap(), cmsRelvalreport::exit, create_public_lumi_plots::exp, fCSC, fDT, DTTopology::firstChannel(), 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, GetRecoTauVFromDQM_MC_cff::next, getHLTPrescaleColumns::path, SteppingHelixStateInfo::path(), PV3DBase< T, PVType, FrameType >::phi(), pi, sysUtil::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(), edm::Event::streamID(), 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().

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

Definition at line 525 of file MuonSimHitProducer.cc.

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

Referenced by MuonSimHitProducer().

527  {
528  // Muons
529  theSimModuleLabel_ = fastMuons.getParameter<std::string>("simModuleLabel");
530  theSimModuleProcess_ = fastMuons.getParameter<std::string>("simModuleProcess");
531  theTrkModuleLabel_ = fastMuons.getParameter<std::string>("trackModuleLabel");
532  std::vector<double> simHitIneffDT = fastMuons.getParameter<std::vector<double> >("simHitDTIneffParameters");
533  std::vector<double> simHitIneffCSC = fastMuons.getParameter<std::vector<double> >("simHitCSCIneffParameters");
534  kDT = simHitIneffDT[0];
535  fDT = simHitIneffDT[1];
536  kCSC = simHitIneffCSC[0];
537  fCSC = simHitIneffCSC[1];
538 
539  // Tracks
540  fullPattern_ = fastTracks.getUntrackedParameter<bool>("FullPatternRecognition");
541 
542  // The following should be on LogInfo
543  // std::cout << " MUON SIM HITS: FastSimulation parameters " << std::endl;
544  // std::cout << " ============================================== " << std::endl;
545  // if ( fullPattern_ )
546  // std::cout << " The FULL pattern recognition option is turned ON" << std::endl;
547  // else
548  // std::cout << " The FAST tracking option is turned ON" << std::endl;
549 
550  // Material Effects
551  theMaterialEffects = 0;
552  if ( matEff.getParameter<bool>("PairProduction") ||
553  matEff.getParameter<bool>("Bremsstrahlung") ||
554  matEff.getParameter<bool>("MuonBremsstrahlung") ||
555  matEff.getParameter<bool>("EnergyLoss") ||
556  matEff.getParameter<bool>("MultipleScattering") )
557  theMaterialEffects = new MaterialEffects(matEff);
558 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
MaterialEffects * theMaterialEffects
std::string theSimModuleProcess_
std::string theSimModuleLabel_
std::string theTrkModuleLabel_

Member Data Documentation

const CSCGeometry* MuonSimHitProducer::cscGeom
private

Definition at line 70 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 69 of file MuonSimHitProducer.h.

Referenced by beginRun(), and produce().

double MuonSimHitProducer::fCSC
private

Definition at line 87 of file MuonSimHitProducer.h.

Referenced by produce(), and readParameters().

double MuonSimHitProducer::fDT
private

Definition at line 85 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 86 of file MuonSimHitProducer.h.

Referenced by produce(), and readParameters().

double MuonSimHitProducer::kDT
private

Definition at line 84 of file MuonSimHitProducer.h.

Referenced by produce(), and readParameters().

const MagneticField* MuonSimHitProducer::magfield
private

Definition at line 68 of file MuonSimHitProducer.h.

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

const Propagator* MuonSimHitProducer::propagatorWithMaterial
private

Definition at line 72 of file MuonSimHitProducer.h.

Referenced by beginRun(), and produce().

Propagator* MuonSimHitProducer::propagatorWithoutMaterial
private

Definition at line 73 of file MuonSimHitProducer.h.

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

const RPCGeometry* MuonSimHitProducer::rpcGeom
private

Definition at line 71 of file MuonSimHitProducer.h.

Referenced by beginRun(), and produce().

Chi2MeasurementEstimator MuonSimHitProducer::theEstimator
private

Definition at line 66 of file MuonSimHitProducer.h.

Referenced by produce().

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

Definition at line 65 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().