CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
fastsim::TrackerSimHitProducer Class Reference

Produces SimHits in the tracker layers. More...

Inheritance diagram for fastsim::TrackerSimHitProducer:
fastsim::InteractionModel

Public Member Functions

std::pair< double, std::unique_ptr< PSimHit > > createHitOnDetector (const TrajectoryStateOnSurface &particle, int pdgId, double layerThickness, double eLoss, int simTrackId, const GeomDet &detector, GlobalPoint &refPos)
 Helper funtion to create the actual SimHit on a detector (sub-) module. More...
 
void interact (Particle &particle, const SimplifiedGeometry &layer, std::vector< std::unique_ptr< Particle > > &secondaries, const RandomEngineAndDistribution &random) override
 Perform the interaction. More...
 
void registerProducts (edm::ProducerBase &producer) const override
 Register the SimHit collection. More...
 
void storeProducts (edm::Event &iEvent) override
 Store the SimHit collection. More...
 
 TrackerSimHitProducer (const std::string &name, const edm::ParameterSet &cfg)
 Constructor. More...
 
 ~TrackerSimHitProducer () override
 Default destructor. More...
 
- Public Member Functions inherited from fastsim::InteractionModel
const std::string getName ()
 Return (unique) name of this interaction. More...
 
 InteractionModel (std::string name)
 Constructor. More...
 
virtual ~InteractionModel ()
 Default destructor. More...
 

Private Attributes

bool doHitsFromInboundParticles_
 If not set, incoming particles (negative speed relative to center of detector) don't create a SimHits since reconstruction anyways not possible. More...
 
double minMomentum_
 Set the minimal momentum of incoming particle. More...
 
const double onSurfaceTolerance_
 Max distance between particle and active (sub-) module. Otherwise particle has to be propagated. More...
 
std::unique_ptr< edm::PSimHitContainersimHitContainer_
 The SimHit. More...
 

Detailed Description

Produces SimHits in the tracker layers.

Also assigns the energy loss of the particle (ionization) with the SimHit. Furthermore, SimHits from different SubModules have to be sorted by their occurance!

See also
EnergyLoss

Definition at line 55 of file TrackerSimHitProducer.cc.

Constructor & Destructor Documentation

fastsim::TrackerSimHitProducer::TrackerSimHitProducer ( const std::string &  name,
const edm::ParameterSet cfg 
)

Constructor.

Definition at line 102 of file TrackerSimHitProducer.cc.

References doHitsFromInboundParticles_, edm::ParameterSet::getParameter(), and minMomentum_.

104  , onSurfaceTolerance_(0.01)
106 {
107  // Set the minimal momentum
108  minMomentum_ = cfg.getParameter<double>("minMomentumCut");
109  // - if not set, particles from outside the beampipe with a negative speed in R direction are propagated but no SimHits
110  // - particle with positive (negative) z and negative (positive) speed in z direction: no SimHits
111  // -> this is not neccesary since a track reconstruction is not possible in this case anyways
112  doHitsFromInboundParticles_ = cfg.getParameter<bool>("doHitsFromInboundParticles");
113 }
T getParameter(std::string const &) const
std::unique_ptr< edm::PSimHitContainer > simHitContainer_
The SimHit.
Base class for any interaction model between a particle and a tracker layer.
bool doHitsFromInboundParticles_
If not set, incoming particles (negative speed relative to center of detector) don&#39;t create a SimHits...
double minMomentum_
Set the minimal momentum of incoming particle.
std::vector< PSimHit > PSimHitContainer
const double onSurfaceTolerance_
Max distance between particle and active (sub-) module. Otherwise particle has to be propagated...
fastsim::TrackerSimHitProducer::~TrackerSimHitProducer ( )
inlineoverride

Member Function Documentation

std::pair< double, std::unique_ptr< PSimHit > > fastsim::TrackerSimHitProducer::createHitOnDetector ( const TrajectoryStateOnSurface particle,
int  pdgId,
double  layerThickness,
double  eLoss,
int  simTrackId,
const GeomDet detector,
GlobalPoint refPos 
)

Helper funtion to create the actual SimHit on a detector (sub-) module.

Parameters
particleRepresentation of the particle's trajectory
pdgIdThe pdgId of the particle
layerThicknessThe thickness of the layer (in radLengths)
eLossThe energy that particle deposited in the detector (lost via ionisation)
simTrackIdThe SimTrack this hit should be assigned to
detectorThe detector element that is hit
refPosReference position that is used to sort the hits if layer consists of sub-detectors \ return Returns the SimHit and the distance relative to refPos since hits have to be ordered (in time) afterwards.

Definition at line 246 of file TrackerSimHitProducer.cc.

References anyDirection, PV3DBase< T, PVType, FrameType >::basicVector(), Surface::bounds(), DEFINE_EDM_PLUGIN, mps_splice::entry, cmsRelvalreport::exit, GeomDet::geographicalId(), TrajectoryStateOnSurface::globalMomentum(), TrajectoryStateOnSurface::globalPosition(), Bounds::length(), TrajectoryStateOnSurface::localMomentum(), TrajectoryStateOnSurface::localPosition(), PV3DBase< T, PVType, FrameType >::mag(), mag(), onSurfaceTolerance_, callgraph::path, common_cff::pdgId, PV3DBase< T, PVType, FrameType >::perp(), PV3DBase< T, PVType, FrameType >::phi(), GloballyPositioned< T >::position(), DetId::rawId(), fastsim::Constants::speedOfLight, GeomDet::surface(), PV3DBase< T, PVType, FrameType >::theta(), Bounds::thickness(), Surface::toGlobal(), GeomDet::toLocal(), TrajectoryStateOnSurface::transverseCurvature(), Vector3DBase< T, FrameTag >::unit(), Bounds::width(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

Referenced by interact(), and ~TrackerSimHitProducer().

247 {
248  // determine position and momentum of particle in the coordinate system of the detector
249  LocalPoint localPosition;
250  LocalVector localMomentum;
251 
252  // if the particle is close enough, no further propagation is needed
253  if(fabs(detector.toLocal(particle.globalPosition()).z()) < onSurfaceTolerance_)
254  {
255  localPosition = particle.localPosition();
256  localMomentum = particle.localMomentum();
257  }
258  // else, propagate
259  else
260  {
261  // find crossing of particle with detector layer
263  particle.globalMomentum().basicVector(),
264  particle.transverseCurvature(),
265  anyDirection);
266  std::pair<bool,double> path = crossing.pathLength(detector.surface());
267  // case propagation succeeds
268  if (path.first)
269  {
270  localPosition = detector.toLocal( GlobalPoint( crossing.position(path.second)));
271  localMomentum = detector.toLocal( GlobalVector( crossing.direction(path.second)));
272  localMomentum = localMomentum.unit() * particle.localMomentum().mag();
273  }
274  // case propagation fails
275  else
276  {
277  return std::pair<double, std::unique_ptr<PSimHit>>(0, std::unique_ptr<PSimHit>(nullptr));
278  }
279  }
280 
281  // find entry and exit point of particle in detector
282  const Plane& detectorPlane = detector.surface();
283  double halfThick = 0.5 * detectorPlane.bounds().thickness();
284  double pZ = localMomentum.z();
285  LocalPoint entry = localPosition + (-halfThick / pZ) * localMomentum;
286  LocalPoint exit = localPosition + (halfThick / pZ) * localMomentum;
287  double tof = particle.globalPosition().mag() / fastsim::Constants::speedOfLight ; // in nanoseconds
288 
289  // make sure the simhit is physically on the module
290  double boundX = detectorPlane.bounds().width() / 2.;
291  double boundY = detectorPlane.bounds().length() / 2.;
292  // Special treatment for TID and TEC trapeziodal modules
293  unsigned subdet = DetId(detector.geographicalId()).subdetId();
294  if (subdet == 4 || subdet == 6)
295  boundX *= 1. - localPosition.y() / detectorPlane.position().perp();
296  if(fabs(localPosition.x()) > boundX || fabs(localPosition.y()) > boundY )
297  {
298  return std::pair<double, std::unique_ptr<PSimHit>>(0, std::unique_ptr<PSimHit>(nullptr));
299  }
300 
301  // Create the hit: the energy loss rescaled to the module thickness
302  // Total thickness is in radiation lengths, 1 radlen = 9.36 cm
303  // Sensitive module thickness is about 30 microns larger than
304  // the module thickness itself
305  eLoss *= (2. * halfThick - 0.003) / (9.36 * layerThickness);
306 
307  // Position of the hit in global coordinates
308  GlobalPoint hitPos(detector.surface().toGlobal(localPosition));
309 
310  return std::pair<double, std::unique_ptr<PSimHit>>((hitPos-refPos).mag(),
311  std::unique_ptr<PSimHit>(new PSimHit(entry, exit, localMomentum.mag(), tof, eLoss, pdgId,
312  detector.geographicalId().rawId(), simTrackId,
313  localMomentum.theta(),
314  localMomentum.phi())));
315 }
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:106
virtual float length() const =0
T perp() const
Definition: PV3DBase.h:72
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
T y() const
Definition: PV3DBase.h:63
GlobalPoint globalPosition() const
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
Definition: GeomDet.h:69
static double constexpr speedOfLight
Speed of light [cm / ns].
Definition: Constants.h:16
const Bounds & bounds() const
Definition: Surface.h:120
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:42
Definition: Plane.h:17
Geom::Theta< T > theta() const
Definition: PV3DBase.h:75
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
LocalVector localMomentum() const
virtual float width() const =0
T mag() const
Definition: PV3DBase.h:67
T z() const
Definition: PV3DBase.h:64
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:79
Vector3DBase unit() const
Definition: Vector3DBase.h:57
Definition: DetId.h:18
virtual float thickness() const =0
GlobalVector globalMomentum() const
T x() const
Definition: PV3DBase.h:62
const BasicVectorType & basicVector() const
Definition: PV3DBase.h:56
const PositionType & position() const
Global3DVector GlobalVector
Definition: GlobalVector.h:10
const double onSurfaceTolerance_
Max distance between particle and active (sub-) module. Otherwise particle has to be propagated...
void fastsim::TrackerSimHitProducer::interact ( Particle particle,
const SimplifiedGeometry layer,
std::vector< std::unique_ptr< Particle > > &  secondaries,
const RandomEngineAndDistribution random 
)
overridevirtual

Perform the interaction.

Parameters
particleThe particle that interacts with the matter.
layerThe detector layer that interacts with the particle.
secondariesParticles that are produced in the interaction (if any).
randomThe Random Engine.

Implements fastsim::InteractionModel.

Definition at line 126 of file TrackerSimHitProducer.cc.

References anyDirection, fastsim::Particle::charge(), GeometricSearchDet::compatibleDets(), GeomDet::components(), createHitOnDetector(), gamEcalExtractorBlocks_cff::detector, doHitsFromInboundParticles_, fastsim::SimplifiedGeometry::getDetLayer(), fastsim::Particle::getEnergyDeposit(), fastsim::SimplifiedGeometry::getMagneticFieldZ(), fastsim::Particle::getMotherSimTrackIndex(), fastsim::SimplifiedGeometry::getThickness(), fastsim::SimplifiedGeometry::isForward(), GeomDet::isLeaf(), fastsim::Particle::isLooper(), seedCreatorFromRegionConsecutiveHitsEDProducer_cff::magneticField, genParticles_cff::map, minMomentum_, fastsim::Particle::momentum(), eostools::move(), common_cff::pdgId, fastsim::Particle::pdgId(), fastsim::Particle::position(), position, PhotonConversionTrajectorySeedProducerFromQuadruplets_cfi::propagator, fastsim::Particle::setEnergyDeposit(), fastsim::Particle::setLooper(), simHitContainer_, fastsim::Particle::simTrackIndex(), GeometricSearchDet::surface(), and Surface::tangentPlane().

Referenced by ~TrackerSimHitProducer().

127 {
128  // the energy deposit in the layer
129  double energyDeposit = particle.getEnergyDeposit();
130  particle.setEnergyDeposit(0);
131 
132  //
133  // don't save hits from particle that did a loop or are inbound (coming from the outer part of the tracker, going towards the center)
134  //
136  if(particle.isLooper()){
137  return;
138  }
139  if(particle.momentum().X()*particle.position().X() + particle.momentum().Y()*particle.position().Y() < 0){
140  particle.setLooper();
141  return;
142  }
143  if(layer.isForward() && ((particle.position().Z() > 0 && particle.momentum().Z() < 0) || (particle.position().Z() < 0 && particle.momentum().Z() > 0))){
144  particle.setLooper();
145  return;
146  }
147  }
148 
149  //
150  // check that layer has tracker modules
151  //
152  if(!layer.getDetLayer())
153  {
154  return;
155  }
156 
157  //
158  // no material
159  //
160  if(layer.getThickness(particle.position(), particle.momentum()) < 1E-10)
161  {
162  return;
163  }
164 
165  //
166  // only charged particles
167  //
168  if(particle.charge()==0)
169  {
170  return;
171  }
172 
173  //
174  // save hit only if momentum higher than threshold
175  //
176  if(particle.momentum().Perp2() < minMomentum_*minMomentum_)
177  {
178  return;
179  }
180 
181  // create the trajectory of the particle
182  UniformMagneticField magneticField(layer.getMagneticFieldZ(particle.position()));
183  GlobalPoint position( particle.position().X(), particle.position().Y(), particle.position().Z());
184  GlobalVector momentum( particle.momentum().Px(), particle.momentum().Py(), particle.momentum().Pz());
185  auto plane = layer.getDetLayer()->surface().tangentPlane(position);
186  TrajectoryStateOnSurface trajectory(GlobalTrajectoryParameters( position, momentum, TrackCharge( particle.charge()), &magneticField), *plane);
187 
188  // find detectors compatible with the particle's trajectory
191  std::vector<DetWithState> compatibleDetectors = layer.getDetLayer()->compatibleDets(trajectory, propagator, est);
192 
194  // You have to sort the simHits in the order they occur!
196 
197  // The old algorithm (sorting by distance to IP) doesn't seem to make sense to me (what if particle moves inwards?)
198 
199  // Detector layers have to be sorted by proximity to particle.position
200  // Propagate particle backwards a bit to make sure it's outside any components
201  std::map<double, std::unique_ptr<PSimHit>> distAndHits;
202  // Position relative to which the hits should be sorted
203  GlobalPoint positionOutside(particle.position().x()-particle.momentum().x()/particle.momentum().P()*10.,
204  particle.position().y()-particle.momentum().y()/particle.momentum().P()*10.,
205  particle.position().z()-particle.momentum().z()/particle.momentum().P()*10.);
206 
207  // FastSim: cheat tracking -> assign hits to closest charged daughter if particle decays
208  int pdgId = particle.pdgId();
209  int simTrackId = particle.getMotherSimTrackIndex() >=0 ? particle.getMotherSimTrackIndex() : particle.simTrackIndex();
210 
211  // loop over the compatible detectors
212  for (const auto & detectorWithState : compatibleDetectors)
213  {
214  const GeomDet & detector = *detectorWithState.first;
215  const TrajectoryStateOnSurface & particleState = detectorWithState.second;
216 
217  // if the detector has no components
218  if(detector.isLeaf())
219  {
220  std::pair<double, std::unique_ptr<PSimHit>> hitPair = createHitOnDetector(particleState, pdgId, layer.getThickness(particle.position()), energyDeposit, simTrackId, detector, positionOutside);
221  if(hitPair.second){
222  distAndHits.insert(distAndHits.end(), std::move(hitPair));
223  }
224  }
225  else
226  {
227  // if the detector has components
228  for(const auto component : detector.components())
229  {
230  std::pair<double, std::unique_ptr<PSimHit>> hitPair = createHitOnDetector(particleState, pdgId,layer.getThickness(particle.position()), energyDeposit, simTrackId, *component, positionOutside);
231  if(hitPair.second){
232  distAndHits.insert(distAndHits.end(), std::move(hitPair));
233  }
234  }
235  }
236  }
237 
238  // Fill simHitContainer
239  for(std::map<double, std::unique_ptr<PSimHit>>::const_iterator it = distAndHits.begin(); it != distAndHits.end(); it++){
240  simHitContainer_->push_back(*(it->second));
241  }
242 
243 }
double Px() const
Definition: Particle.h:57
double Py() const
Definition: Particle.h:59
std::unique_ptr< edm::PSimHitContainer > simHitContainer_
The SimHit.
double X() const
Definition: Particle.h:49
double Z() const
Definition: Particle.h:53
int TrackCharge
Definition: TrackCharge.h:4
bool doHitsFromInboundParticles_
If not set, incoming particles (negative speed relative to center of detector) don&#39;t create a SimHits...
std::pair< double, std::unique_ptr< PSimHit > > createHitOnDetector(const TrajectoryStateOnSurface &particle, int pdgId, double layerThickness, double eLoss, int simTrackId, const GeomDet &detector, GlobalPoint &refPos)
Helper funtion to create the actual SimHit on a detector (sub-) module.
double Pz() const
Definition: Particle.h:61
virtual bool isLeaf() const
is a Unit
Definition: GeomDet.h:85
virtual std::vector< const GeomDet * > components() const
Returns direct components, if any.
Definition: GeomDet.h:88
double minMomentum_
Set the minimal momentum of incoming particle.
static int position[264][3]
Definition: ReadPGInfo.cc:509
def move(src, dest)
Definition: eostools.py:510
double Y() const
Definition: Particle.h:51
void fastsim::TrackerSimHitProducer::registerProducts ( edm::ProducerBase producer) const
overridevirtual

Register the SimHit collection.

Reimplemented from fastsim::InteractionModel.

Definition at line 115 of file TrackerSimHitProducer.cc.

References edm::ProductRegistryHelper::produces().

Referenced by ~TrackerSimHitProducer().

116 {
117  producer.produces<edm::PSimHitContainer>("TrackerHits");
118 }
BranchAliasSetterT< ProductType > produces()
declare what type of product will make and with which optional label
std::vector< PSimHit > PSimHitContainer
void fastsim::TrackerSimHitProducer::storeProducts ( edm::Event iEvent)
overridevirtual

Store the SimHit collection.

Reimplemented from fastsim::InteractionModel.

Definition at line 120 of file TrackerSimHitProducer.cc.

References eostools::move(), edm::Event::put(), and simHitContainer_.

Referenced by ~TrackerSimHitProducer().

121 {
122  iEvent.put(std::move(simHitContainer_), "TrackerHits");
124 }
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:127
std::unique_ptr< edm::PSimHitContainer > simHitContainer_
The SimHit.
std::vector< PSimHit > PSimHitContainer
def move(src, dest)
Definition: eostools.py:510

Member Data Documentation

bool fastsim::TrackerSimHitProducer::doHitsFromInboundParticles_
private

If not set, incoming particles (negative speed relative to center of detector) don't create a SimHits since reconstruction anyways not possible.

Definition at line 96 of file TrackerSimHitProducer.cc.

Referenced by interact(), and TrackerSimHitProducer().

double fastsim::TrackerSimHitProducer::minMomentum_
private

Set the minimal momentum of incoming particle.

Definition at line 95 of file TrackerSimHitProducer.cc.

Referenced by interact(), and TrackerSimHitProducer().

const double fastsim::TrackerSimHitProducer::onSurfaceTolerance_
private

Max distance between particle and active (sub-) module. Otherwise particle has to be propagated.

Definition at line 93 of file TrackerSimHitProducer.cc.

Referenced by createHitOnDetector().

std::unique_ptr<edm::PSimHitContainer> fastsim::TrackerSimHitProducer::simHitContainer_
private

The SimHit.

Definition at line 94 of file TrackerSimHitProducer.cc.

Referenced by interact(), and storeProducts().